]> https://gitweb.dealii.org/ - dealii.git/commitdiff
Early deprecate functions in DoFTools 13365/head
authorPeter Munch <peterrmuench@gmail.com>
Sat, 12 Feb 2022 06:45:51 +0000 (07:45 +0100)
committerPeter Munch <peterrmuench@gmail.com>
Sun, 13 Feb 2022 13:51:53 +0000 (14:51 +0100)
25 files changed:
doc/news/changes/incompatibilities/20220212Munch [new file with mode: 0644]
examples/step-16/step-16.cc
examples/step-18/step-18.cc
examples/step-32/step-32.cc
examples/step-37/doc/results.dox
examples/step-37/step-37.cc
examples/step-40/step-40.cc
examples/step-42/step-42.cc
examples/step-45/step-45.cc
examples/step-48/step-48.cc
examples/step-50/step-50.cc
examples/step-55/step-55.cc
examples/step-62/step-62.cc
examples/step-63/step-63.cc
examples/step-64/step-64.cu
examples/step-66/step-66.cc
examples/step-68/step-68.cc
examples/step-69/step-69.cc
examples/step-70/step-70.cc
examples/step-75/step-75.cc
include/deal.II/dofs/dof_tools.h
include/deal.II/matrix_free/cuda_matrix_free.templates.h
include/deal.II/multigrid/mg_constrained_dofs.h
source/dofs/dof_tools.cc
source/dofs/dof_tools.inst.in

diff --git a/doc/news/changes/incompatibilities/20220212Munch b/doc/news/changes/incompatibilities/20220212Munch
new file mode 100644 (file)
index 0000000..ed55891
--- /dev/null
@@ -0,0 +1,9 @@
+Deprecated: We have introduced new versions of
+DoFTools::extract_locally_active_dofs(),
+DoFTools::extract_locally_active_level_dofs(),
+DoFTools::extract_locally_relevant_dofs(), and
+DoFTools::extract_locally_relevant_level_dofs().
+These versions return the index sets directly.
+The old versions have been early deprecated.
+<br>
+(Peter Munch, 2022/02/23)
index cd58b79527f90a722004a1c62b0e75e5fc635a09..ea76edf5ef218e4ee02dbc25fe3fd0eb905272df 100644 (file)
@@ -432,10 +432,8 @@ namespace Step16
     std::vector<AffineConstraints<double>> boundary_constraints(n_levels);
     for (unsigned int level = 0; level < n_levels; ++level)
       {
-        IndexSet dofset;
-        DoFTools::extract_locally_relevant_level_dofs(dof_handler,
-                                                      level,
-                                                      dofset);
+        const IndexSet dofset =
+          DoFTools::extract_locally_relevant_level_dofs(dof_handler, level);
         boundary_constraints[level].reinit(dofset);
         boundary_constraints[level].add_lines(
           mg_constrained_dofs.get_refinement_edge_indices(level));
index ba907366d17329cad3bdd7cc29838da54fcc2a7f..59d96acd82b8afa10b4f2deeac06f172c502f060 100644 (file)
@@ -800,7 +800,8 @@ namespace Step18
   {
     dof_handler.distribute_dofs(fe);
     locally_owned_dofs = dof_handler.locally_owned_dofs();
-    DoFTools::extract_locally_relevant_dofs(dof_handler, locally_relevant_dofs);
+    locally_relevant_dofs =
+      DoFTools::extract_locally_relevant_dofs(dof_handler);
 
     // The next step is to set up constraints due to hanging nodes. This has
     // been handled many times before:
index 496ffe4ffd19ec72e6f4246b6843e2c3f3f7d2c4..cfb968f659421819e26555240a775077652ff45c 100644 (file)
@@ -1864,16 +1864,16 @@ namespace Step32
       stokes_partitioning.push_back(stokes_index_set.get_view(0, n_u));
       stokes_partitioning.push_back(stokes_index_set.get_view(n_u, n_u + n_p));
 
-      DoFTools::extract_locally_relevant_dofs(stokes_dof_handler,
-                                              stokes_relevant_set);
+      stokes_relevant_set =
+        DoFTools::extract_locally_relevant_dofs(stokes_dof_handler);
       stokes_relevant_partitioning.push_back(
         stokes_relevant_set.get_view(0, n_u));
       stokes_relevant_partitioning.push_back(
         stokes_relevant_set.get_view(n_u, n_u + n_p));
 
       temperature_partitioning = temperature_dof_handler.locally_owned_dofs();
-      DoFTools::extract_locally_relevant_dofs(
-        temperature_dof_handler, temperature_relevant_partitioning);
+      temperature_relevant_partitioning =
+        DoFTools::extract_locally_relevant_dofs(temperature_dof_handler);
     }
 
     // Following this, we can compute constraints for the solution vectors,
index e02ffd26c1ba1dba478e3f1f2512686d16cb0259..7ef0edf09bbfe77337c6c7040a985b0c02f3f10e 100644 (file)
@@ -427,8 +427,7 @@ temporary vector and LinearAlgebra::distributed::Vector::copy_locally_owned_data
 as shown below.
 
 @code
-IndexSet locally_relevant_dofs;
-DoFTools::extract_locally_relevant_dofs(dof_handler, locally_relevant_dofs);
+const IndexSet locally_relevant_dofs = DoFTools::extract_locally_relevant_dofs(dof_handler);
 LinearAlgebra::distributed::Vector<double> copy_vec(solution);
 solution.reinit(dof_handler.locally_owned_dofs(),
                 locally_relevant_dofs,
index 134fca188060ebaedac6cb6cd7fc838fc838a620..69f01c762a48a89676157fbf026f4b534832ee92 100644 (file)
@@ -792,8 +792,8 @@ namespace Step37
     pcout << "Number of degrees of freedom: " << dof_handler.n_dofs()
           << std::endl;
 
-    IndexSet locally_relevant_dofs;
-    DoFTools::extract_locally_relevant_dofs(dof_handler, locally_relevant_dofs);
+    const IndexSet locally_relevant_dofs =
+      DoFTools::extract_locally_relevant_dofs(dof_handler);
 
     constraints.clear();
     constraints.reinit(locally_relevant_dofs);
@@ -852,10 +852,8 @@ namespace Step37
 
     for (unsigned int level = 0; level < nlevels; ++level)
       {
-        IndexSet relevant_dofs;
-        DoFTools::extract_locally_relevant_level_dofs(dof_handler,
-                                                      level,
-                                                      relevant_dofs);
+        const IndexSet relevant_dofs =
+          DoFTools::extract_locally_relevant_level_dofs(dof_handler, level);
         AffineConstraints<double> level_constraints;
         level_constraints.reinit(relevant_dofs);
         level_constraints.add_lines(
index 074e9c11e8501d97968f0da46aaba97066395dc4..1f61a65f72dbbed4e4b687b628298915a286062f 100644 (file)
@@ -256,7 +256,8 @@ namespace Step40
     // around the locally owned cells; we need all of these degrees of
     // freedom, for example, to estimate the error on the local cells).
     locally_owned_dofs = dof_handler.locally_owned_dofs();
-    DoFTools::extract_locally_relevant_dofs(dof_handler, locally_relevant_dofs);
+    locally_relevant_dofs =
+      DoFTools::extract_locally_relevant_dofs(dof_handler);
 
     // Next, let us initialize the solution and right hand side vectors. As
     // mentioned above, the solution vector we seek does not only store
index 5db3a410baf36bac191bad4d7294c21b447dd064..cfa0dd7822c6df291989129257102b9caf28b5c7 100644 (file)
@@ -973,9 +973,8 @@ namespace Step42
       dof_handler.distribute_dofs(fe);
 
       locally_owned_dofs = dof_handler.locally_owned_dofs();
-      locally_relevant_dofs.clear();
-      DoFTools::extract_locally_relevant_dofs(dof_handler,
-                                              locally_relevant_dofs);
+      locally_relevant_dofs =
+        DoFTools::extract_locally_relevant_dofs(dof_handler);
     }
 
     /* setup hanging nodes and Dirichlet constraints */
index 488b526bc97f0683666b32515e33dd200ddf64f1..3a3cc6f1af3b942a512d272c05fbba56283ebd4a 100644 (file)
@@ -374,9 +374,8 @@ namespace Step45
       owned_partitioning.push_back(locally_owned_dofs.get_view(n_u, n_u + n_p));
 
       relevant_partitioning.clear();
-      IndexSet locally_relevant_dofs;
-      DoFTools::extract_locally_relevant_dofs(dof_handler,
-                                              locally_relevant_dofs);
+      const IndexSet locally_relevant_dofs =
+        DoFTools::extract_locally_relevant_dofs(dof_handler);
       relevant_partitioning.push_back(locally_relevant_dofs.get_view(0, n_u));
       relevant_partitioning.push_back(
         locally_relevant_dofs.get_view(n_u, n_u + n_p));
index 8daa6d4cc6368307a418aacecae7272aab9fc68c..1255f5f61f83bd5ca1668dacf4d5012a131a1525 100644 (file)
@@ -427,7 +427,8 @@ namespace Step48
     // access in MPI-local numbers that need to match between the vector and
     // MatrixFree), so we just ask it to initialize the vectors to be sure the
     // ghost exchange is properly handled.
-    DoFTools::extract_locally_relevant_dofs(dof_handler, locally_relevant_dofs);
+    locally_relevant_dofs =
+      DoFTools::extract_locally_relevant_dofs(dof_handler);
     constraints.clear();
     constraints.reinit(locally_relevant_dofs);
     DoFTools::make_hanging_node_constraints(dof_handler, constraints);
index 413808475455aa99414aad22b04df8d2f8f64259..e118e32800f909c70e8af23b47a5f0581748fc38 100644 (file)
@@ -487,8 +487,8 @@ void LaplaceProblem<dim, degree>::setup_system()
 
   dof_handler.distribute_dofs(fe);
 
-  DoFTools::extract_locally_relevant_dofs(dof_handler, locally_relevant_dofs);
-  locally_owned_dofs = dof_handler.locally_owned_dofs();
+  locally_relevant_dofs = DoFTools::extract_locally_relevant_dofs(dof_handler);
+  locally_owned_dofs    = dof_handler.locally_owned_dofs();
 
   solution.reinit(locally_owned_dofs, mpi_communicator);
   right_hand_side.reinit(locally_owned_dofs, mpi_communicator);
@@ -594,10 +594,9 @@ void LaplaceProblem<dim, degree>::setup_multigrid()
 
           for (unsigned int level = 0; level < n_levels; ++level)
             {
-              IndexSet relevant_dofs;
-              DoFTools::extract_locally_relevant_level_dofs(dof_handler,
-                                                            level,
-                                                            relevant_dofs);
+              const IndexSet relevant_dofs =
+                DoFTools::extract_locally_relevant_level_dofs(dof_handler,
+                                                              level);
               AffineConstraints<double> level_constraints;
               level_constraints.reinit(relevant_dofs);
               level_constraints.add_lines(
@@ -642,10 +641,9 @@ void LaplaceProblem<dim, degree>::setup_multigrid()
 
           for (unsigned int level = 0; level < n_levels; ++level)
             {
-              IndexSet dof_set;
-              DoFTools::extract_locally_relevant_level_dofs(dof_handler,
-                                                            level,
-                                                            dof_set);
+              const IndexSet dof_set =
+                DoFTools::extract_locally_relevant_level_dofs(dof_handler,
+                                                              level);
 
               {
 #ifdef USE_PETSC_LA
@@ -828,10 +826,8 @@ void LaplaceProblem<dim, degree>::assemble_multigrid()
     triangulation.n_global_levels());
   for (unsigned int level = 0; level < triangulation.n_global_levels(); ++level)
     {
-      IndexSet dof_set;
-      DoFTools::extract_locally_relevant_level_dofs(dof_handler,
-                                                    level,
-                                                    dof_set);
+      const IndexSet dof_set =
+        DoFTools::extract_locally_relevant_level_dofs(dof_handler, level);
       boundary_constraints[level].reinit(dof_set);
       boundary_constraints[level].add_lines(
         mg_constrained_dofs.get_refinement_edge_indices(level));
index b5d9256372856087f1efa3faa042afaab50a1487..854342e0268f86ac40e72dbc9a32a3a8c61f88ec 100644 (file)
@@ -382,8 +382,8 @@ namespace Step55
     owned_partitioning[1] =
       dof_handler.locally_owned_dofs().get_view(n_u, n_u + n_p);
 
-    IndexSet locally_relevant_dofs;
-    DoFTools::extract_locally_relevant_dofs(dof_handler, locally_relevant_dofs);
+    const IndexSet locally_relevant_dofs =
+      DoFTools::extract_locally_relevant_dofs(dof_handler);
     relevant_partitioning.resize(2);
     relevant_partitioning[0] = locally_relevant_dofs.get_view(0, n_u);
     relevant_partitioning[1] = locally_relevant_dofs.get_view(n_u, n_u + n_p);
index 9e2614eec6bba532421abb1b67ff61e85895da4f..1d8a9572caf9de4db70ba1a10471bc55e1faf4f1 100644 (file)
@@ -748,7 +748,8 @@ namespace step62
     dof_handler.distribute_dofs(fe);
 
     locally_owned_dofs = dof_handler.locally_owned_dofs();
-    DoFTools::extract_locally_relevant_dofs(dof_handler, locally_relevant_dofs);
+    locally_relevant_dofs =
+      DoFTools::extract_locally_relevant_dofs(dof_handler);
 
     locally_relevant_solution.reinit(locally_owned_dofs,
                                      locally_relevant_dofs,
index 01db361eafc6e13bc9cffd55000a6cf761427405..6ac169c569ced168f7dcc5cae7cd02f8fc87c98e 100644 (file)
@@ -814,9 +814,8 @@ namespace Step63
     for (unsigned int level = 0; level < triangulation.n_global_levels();
          ++level)
       {
-        IndexSet locally_owned_level_dof_indices;
-        DoFTools::extract_locally_relevant_level_dofs(
-          dof_handler, level, locally_owned_level_dof_indices);
+        const IndexSet locally_owned_level_dof_indices =
+          DoFTools::extract_locally_relevant_level_dofs(dof_handler, level);
         boundary_constraints[level].reinit(locally_owned_level_dof_indices);
         boundary_constraints[level].add_lines(
           mg_constrained_dofs.get_refinement_edge_indices(level));
index cb76fea302e9627241e1756955993bd932c65700..bc2c550a3280ede785eb0b54455adea098f373fc 100644 (file)
@@ -397,7 +397,8 @@ namespace Step64
     dof_handler.distribute_dofs(fe);
 
     locally_owned_dofs = dof_handler.locally_owned_dofs();
-    DoFTools::extract_locally_relevant_dofs(dof_handler, locally_relevant_dofs);
+    locally_relevant_dofs =
+      DoFTools::extract_locally_relevant_dofs(dof_handler);
     system_rhs_dev.reinit(locally_owned_dofs, mpi_communicator);
 
     constraints.clear();
index 155f7253ca1a75a7e5af6abcaa30c3d02d392e6b..67fcdc6726610da48e6b0b4242ddd589f26502e6 100644 (file)
@@ -544,8 +544,8 @@ namespace Step66
     dof_handler.distribute_dofs(fe);
     dof_handler.distribute_mg_dofs();
 
-    IndexSet locally_relevant_dofs;
-    DoFTools::extract_locally_relevant_dofs(dof_handler, locally_relevant_dofs);
+    const IndexSet locally_relevant_dofs =
+      DoFTools::extract_locally_relevant_dofs(dof_handler);
 
     constraints.clear();
     constraints.reinit(locally_relevant_dofs);
@@ -593,10 +593,8 @@ namespace Step66
 
     for (unsigned int level = 0; level < nlevels; ++level)
       {
-        IndexSet relevant_dofs;
-        DoFTools::extract_locally_relevant_level_dofs(dof_handler,
-                                                      level,
-                                                      relevant_dofs);
+        const IndexSet relevant_dofs =
+          DoFTools::extract_locally_relevant_level_dofs(dof_handler, level);
 
         AffineConstraints<double> level_constraints;
         level_constraints.reinit(relevant_dofs);
index 1dd1d6d05a9b7451a252c08da89b3e12414526f2..2d7bacd367087b17224a1bbd75ad2436d1110202 100644 (file)
@@ -474,8 +474,8 @@ namespace Step68
   {
     fluid_dh.distribute_dofs(fluid_fe);
     const IndexSet locally_owned_dofs = fluid_dh.locally_owned_dofs();
-    IndexSet       locally_relevant_dofs;
-    DoFTools::extract_locally_relevant_dofs(fluid_dh, locally_relevant_dofs);
+    const IndexSet locally_relevant_dofs =
+      DoFTools::extract_locally_relevant_dofs(fluid_dh);
 
     velocity_field.reinit(locally_owned_dofs,
                           locally_relevant_dofs,
index 29683071bfdb0b5a9cd3115c6be747cf9e7163c4..c4928e5551655d3a3df76cba8c5f603f1ef5a429 100644 (file)
@@ -713,7 +713,7 @@ namespace Step69
       locally_owned   = dof_handler.locally_owned_dofs();
       n_locally_owned = locally_owned.n_elements();
 
-      DoFTools::extract_locally_relevant_dofs(dof_handler, locally_relevant);
+      locally_relevant   = DoFTools::extract_locally_relevant_dofs(dof_handler);
       n_locally_relevant = locally_relevant.n_elements();
 
       partitioner =
index 1833fc9618772c437f254fa33bf3cceed826ad91..8d498123cc22e42a49a0d777df6e07887bfb9025 100644 (file)
@@ -1231,8 +1231,8 @@ namespace Step70
     fluid_owned_dofs[1] =
       fluid_dh.locally_owned_dofs().get_view(n_u, n_u + n_p);
 
-    IndexSet locally_relevant_dofs;
-    DoFTools::extract_locally_relevant_dofs(fluid_dh, locally_relevant_dofs);
+    const IndexSet locally_relevant_dofs =
+      DoFTools::extract_locally_relevant_dofs(fluid_dh);
     fluid_relevant_dofs.resize(2);
     fluid_relevant_dofs[0] = locally_relevant_dofs.get_view(0, n_u);
     fluid_relevant_dofs[1] = locally_relevant_dofs.get_view(n_u, n_u + n_p);
index 1ed9c2d60a2150f6f247d9e51d99a8b8dee46464..94ec5f0a7542af86827f5165ff05793ad6e48d95 100644 (file)
@@ -324,9 +324,8 @@ namespace Step75
     {
       AffineConstraints<number> constraints_without_dbc;
 
-      IndexSet locally_relevant_dofs;
-      DoFTools::extract_locally_relevant_dofs(dof_handler,
-                                              locally_relevant_dofs);
+      const IndexSet locally_relevant_dofs =
+        DoFTools::extract_locally_relevant_dofs(dof_handler);
       constraints_without_dbc.reinit(locally_relevant_dofs);
 
       DoFTools::make_hanging_node_constraints(dof_handler,
@@ -777,9 +776,8 @@ namespace Step75
         const auto &dof_handler = dof_handlers[level];
         auto &      constraint  = constraints[level];
 
-        IndexSet locally_relevant_dofs;
-        DoFTools::extract_locally_relevant_dofs(dof_handler,
-                                                locally_relevant_dofs);
+        const IndexSet locally_relevant_dofs =
+          DoFTools::extract_locally_relevant_dofs(dof_handler);
         constraint.reinit(locally_relevant_dofs);
 
         DoFTools::make_hanging_node_constraints(dof_handler, constraint);
@@ -1105,7 +1103,8 @@ namespace Step75
     dof_handler.distribute_dofs(fe_collection);
 
     locally_owned_dofs = dof_handler.locally_owned_dofs();
-    DoFTools::extract_locally_relevant_dofs(dof_handler, locally_relevant_dofs);
+    locally_relevant_dofs =
+      DoFTools::extract_locally_relevant_dofs(dof_handler);
 
     locally_relevant_solution.reinit(locally_owned_dofs,
                                      locally_relevant_dofs,
index ac7fd4621bacc63735447405e224493510d41d1d..5372aa365fafbc38a19cbbf92954e986bb8c6354 100644 (file)
@@ -1554,6 +1554,26 @@ namespace DoFTools
    * with the locally owned subdomain id.
    */
   template <int dim, int spacedim>
+  IndexSet
+  extract_locally_active_dofs(const DoFHandler<dim, spacedim> &dof_handler);
+
+  /**
+   * Extract the set of global DoF indices that are active on the current
+   * DoFHandler. For regular DoFHandlers, these are all DoF indices, but for
+   * DoFHandler objects built on parallel::distributed::Triangulation this set
+   * is a superset of DoFHandler::locally_owned_dofs() and contains all DoF
+   * indices that live on all locally owned cells (including on the interface
+   * to ghost cells). However, it does not contain the DoF indices that are
+   * exclusively defined on ghost or artificial cells (see
+   * @ref GlossArtificialCell "the glossary").
+   *
+   * The degrees of freedom identified by this function equal those obtained
+   * from the dof_indices_with_subdomain_association() function when called
+   * with the locally owned subdomain id.
+   *
+   * @deprecated Use the previous function instead.
+   */
+  template <int dim, int spacedim>
   void
   extract_locally_active_dofs(const DoFHandler<dim, spacedim> &dof_handler,
                               IndexSet &                       dof_set);
@@ -1565,6 +1585,20 @@ namespace DoFTools
    * given level.
    */
   template <int dim, int spacedim>
+  IndexSet
+  extract_locally_active_level_dofs(
+    const DoFHandler<dim, spacedim> &dof_handler,
+    const unsigned int               level);
+
+  /**
+   * Same function as above but for a certain (multigrid-)level.
+   * This function returns all DoF indices that live on
+   * all locally owned cells (including on the interface to ghost cells) on the
+   * given level.
+   *
+   * @deprecated Use the previous function instead.
+   */
+  template <int dim, int spacedim>
   void
   extract_locally_active_level_dofs(
     const DoFHandler<dim, spacedim> &dof_handler,
@@ -1581,11 +1615,25 @@ namespace DoFTools
    * @ref GlossArtificialCell "the glossary").
    */
   template <int dim, int spacedim>
+  IndexSet
+  extract_locally_relevant_dofs(const DoFHandler<dim, spacedim> &dof_handler);
+
+  /**
+   * Extract the set of global DoF indices that are active on the current
+   * DoFHandler. For regular DoFHandlers, these are all DoF indices, but for
+   * DoFHandler objects built on parallel::distributed::Triangulation this set
+   * is the union of DoFHandler::locally_owned_dofs() and the DoF indices on
+   * all ghost cells. In essence, it is the DoF indices on all cells that are
+   * not artificial (see
+   * @ref GlossArtificialCell "the glossary").
+   *
+   * @deprecated Use the previous function instead.
+   */
+  template <int dim, int spacedim>
   void
   extract_locally_relevant_dofs(const DoFHandler<dim, spacedim> &dof_handler,
                                 IndexSet &                       dof_set);
 
-
   /**
    * Extract the set of locally owned DoF indices for each component within the
    * mask that are owned by the current  processor. For components disabled by
@@ -1638,10 +1686,21 @@ namespace DoFTools
   locally_relevant_dofs_per_subdomain(
     const DoFHandler<dim, spacedim> &dof_handler);
 
+  /**
+   * Same as extract_locally_relevant_dofs() but for multigrid DoFs for the
+   * given @p level.
+   */
+  template <int dim, int spacedim>
+  IndexSet
+  extract_locally_relevant_level_dofs(
+    const DoFHandler<dim, spacedim> &dof_handler,
+    const unsigned int               level);
 
   /**
    * Same as extract_locally_relevant_dofs() but for multigrid DoFs for the
    * given @p level.
+   *
+   * @deprecated Use the previous function instead.
    */
   template <int dim, int spacedim>
   void
index 08013ef6ba760541f7c2c3bf3c07a5336df856e6..55447ec3311fc83f0a219abb6bd669fa747d1e10 100644 (file)
@@ -1077,8 +1077,8 @@ namespace CUDAWrappers
     IndexSet locally_relevant_dofs;
     if (comm)
       {
-        DoFTools::extract_locally_relevant_dofs(*dof_handler,
-                                                locally_relevant_dofs);
+        locally_relevant_dofs =
+          DoFTools::extract_locally_relevant_dofs(*dof_handler);
         partitioner = std::make_shared<Utilities::MPI::Partitioner>(
           dof_handler->locally_owned_dofs(), locally_relevant_dofs, *comm);
       }
index 6c2a667310ca0d57770b9889751e385e8338c6ed..9d9efbd31767aad3959081028f4b1f4d85aee30d 100644 (file)
@@ -262,8 +262,8 @@ MGConstrainedDoFs::initialize(
         }
       else
         {
-          IndexSet relevant_dofs;
-          DoFTools::extract_locally_relevant_level_dofs(dof, l, relevant_dofs);
+          const IndexSet relevant_dofs =
+            DoFTools::extract_locally_relevant_level_dofs(dof, l);
           level_constraints[l].reinit(relevant_dofs);
         }
 
index 1acf4d5602d531d7a50a68f5b1237a5e42dff9bd..d2dd3e04485179e28ba1071e3198f30d16eec961 100644 (file)
@@ -1048,12 +1048,11 @@ namespace DoFTools
 
 
   template <int dim, int spacedim>
-  void
-  extract_locally_active_dofs(const DoFHandler<dim, spacedim> &dof_handler,
-                              IndexSet &                       dof_set)
+  IndexSet
+  extract_locally_active_dofs(const DoFHandler<dim, spacedim> &dof_handler)
   {
     // collect all the locally owned dofs
-    dof_set = dof_handler.locally_owned_dofs();
+    IndexSet dof_set = dof_handler.locally_owned_dofs();
 
     // add the DoF on the adjacent ghost cells to the IndexSet, cache them
     // in a set. need to check each dof manually because we can't be sure
@@ -1075,19 +1074,30 @@ namespace DoFTools
     dof_set.add_indices(global_dof_indices.begin(), global_dof_indices.end());
 
     dof_set.compress();
+
+    return dof_set;
   }
 
 
 
   template <int dim, int spacedim>
   void
+  extract_locally_active_dofs(const DoFHandler<dim, spacedim> &dof_handler,
+                              IndexSet &                       dof_set)
+  {
+    dof_set = extract_locally_active_dofs(dof_handler);
+  }
+
+
+
+  template <int dim, int spacedim>
+  IndexSet
   extract_locally_active_level_dofs(
     const DoFHandler<dim, spacedim> &dof_handler,
-    IndexSet &                       dof_set,
     const unsigned int               level)
   {
     // collect all the locally owned dofs
-    dof_set = dof_handler.locally_owned_mg_dofs(level);
+    IndexSet dof_set = dof_handler.locally_owned_mg_dofs(level);
 
     // add the DoF on the adjacent ghost cells to the IndexSet, cache them
     // in a set. need to check each dof manually because we can't be sure
@@ -1111,17 +1121,30 @@ namespace DoFTools
     dof_set.add_indices(global_dof_indices.begin(), global_dof_indices.end());
 
     dof_set.compress();
+
+    return dof_set;
   }
 
 
 
   template <int dim, int spacedim>
   void
-  extract_locally_relevant_dofs(const DoFHandler<dim, spacedim> &dof_handler,
-                                IndexSet &                       dof_set)
+  extract_locally_active_level_dofs(
+    const DoFHandler<dim, spacedim> &dof_handler,
+    IndexSet &                       dof_set,
+    const unsigned int               level)
+  {
+    dof_set = extract_locally_active_level_dofs(dof_handler, level);
+  }
+
+
+
+  template <int dim, int spacedim>
+  IndexSet
+  extract_locally_relevant_dofs(const DoFHandler<dim, spacedim> &dof_handler)
   {
     // collect all the locally owned dofs
-    dof_set = dof_handler.locally_owned_dofs();
+    IndexSet dof_set = dof_handler.locally_owned_dofs();
 
     // now add the DoF on the adjacent ghost cells to the IndexSet
 
@@ -1150,19 +1173,30 @@ namespace DoFTools
                         std::unique(dofs_on_ghosts.begin(),
                                     dofs_on_ghosts.end()));
     dof_set.compress();
+
+    return dof_set;
   }
 
 
 
   template <int dim, int spacedim>
   void
+  extract_locally_relevant_dofs(const DoFHandler<dim, spacedim> &dof_handler,
+                                IndexSet &                       dof_set)
+  {
+    dof_set = extract_locally_relevant_dofs(dof_handler);
+  }
+
+
+
+  template <int dim, int spacedim>
+  IndexSet
   extract_locally_relevant_level_dofs(
     const DoFHandler<dim, spacedim> &dof_handler,
-    const unsigned int               level,
-    IndexSet &                       dof_set)
+    const unsigned int               level)
   {
     // collect all the locally owned dofs
-    dof_set = dof_handler.locally_owned_mg_dofs(level);
+    IndexSet dof_set = dof_handler.locally_owned_mg_dofs(level);
 
     // add the DoF on the adjacent ghost cells to the IndexSet
 
@@ -1198,6 +1232,20 @@ namespace DoFTools
                                     dofs_on_ghosts.end()));
 
     dof_set.compress();
+
+    return dof_set;
+  }
+
+
+
+  template <int dim, int spacedim>
+  void
+  extract_locally_relevant_level_dofs(
+    const DoFHandler<dim, spacedim> &dof_handler,
+    const unsigned int               level,
+    IndexSet &                       dof_set)
+  {
+    dof_set = extract_locally_relevant_level_dofs(dof_handler, level);
   }
 
 
index 32898e84bfdef666ddfd5c037c7ead94d331b5f7..c8ae1b9990b33a22988767b1cb3dd3d8b067f3a2 100644 (file)
@@ -19,12 +19,24 @@ for (deal_II_dimension, deal_II_space_dimension : DIMENSIONS)
 #if deal_II_dimension <= deal_II_space_dimension
     namespace DoFTools
     \{
+      template IndexSet
+      extract_locally_relevant_dofs<deal_II_dimension, deal_II_space_dimension>(
+        const DoFHandler<deal_II_dimension, deal_II_space_dimension>
+          &dof_handler);
+
       template void
       extract_locally_relevant_dofs<deal_II_dimension, deal_II_space_dimension>(
         const DoFHandler<deal_II_dimension, deal_II_space_dimension>
           &       dof_handler,
         IndexSet &dof_set);
 
+      template IndexSet
+      extract_locally_relevant_level_dofs<deal_II_dimension,
+                                          deal_II_space_dimension>(
+        const DoFHandler<deal_II_dimension, deal_II_space_dimension>
+          &                dof_handler,
+        const unsigned int level);
+
       template void
       extract_locally_relevant_level_dofs<deal_II_dimension,
                                           deal_II_space_dimension>(
@@ -33,12 +45,24 @@ for (deal_II_dimension, deal_II_space_dimension : DIMENSIONS)
         const unsigned int level,
         IndexSet &         dof_set);
 
+      template IndexSet
+      extract_locally_active_dofs<deal_II_dimension, deal_II_space_dimension>(
+        const DoFHandler<deal_II_dimension, deal_II_space_dimension>
+          &dof_handler);
+
       template void
       extract_locally_active_dofs<deal_II_dimension, deal_II_space_dimension>(
         const DoFHandler<deal_II_dimension, deal_II_space_dimension>
           &       dof_handler,
         IndexSet &dof_set);
 
+      template IndexSet
+      extract_locally_active_level_dofs<deal_II_dimension,
+                                        deal_II_space_dimension>(
+        const DoFHandler<deal_II_dimension, deal_II_space_dimension>
+          &                dof_handler,
+        const unsigned int level);
+
       template void
       extract_locally_active_level_dofs<deal_II_dimension,
                                         deal_II_space_dimension>(

In the beginning the Universe was created. This has made a lot of people very angry and has been widely regarded as a bad move.

Douglas Adams


Typeset in Trocchi and Trocchi Bold Sans Serif.