]> https://gitweb.dealii.org/ - dealii.git/commitdiff
Simplify the access to MG dofs in DoFAccessor
authorMartin Kronbichler <martin.kronbichler@uni-a.de>
Wed, 25 May 2022 22:20:56 +0000 (00:20 +0200)
committerMartin Kronbichler <martin.kronbichler@uni-a.de>
Fri, 3 Jun 2022 15:26:02 +0000 (17:26 +0200)
include/deal.II/dofs/dof_accessor.templates.h
include/deal.II/dofs/dof_handler.h
include/deal.II/dofs/dof_objects.h
source/dofs/dof_handler.cc
source/dofs/dof_handler.inst.in
source/dofs/dof_handler_policy.cc
source/dofs/dof_objects.cc
source/dofs/dof_objects.inst.in

index 9c0203030e3d5785538a112f783c21219a3b6c97..66129d1792a67a8887474a5fd5922c9fb22401a1 100644 (file)
@@ -490,6 +490,48 @@ namespace internal
         return {0, 0};
       }
 
+
+
+      // The next few internal helper functions are needed to support various
+      // DoFIndicesType kinds, e.g. actual vectors of DoFIndices or dummy
+      // types that we use when we only want to work on the internally stored
+      // DoFs and never extract any number.
+      static unsigned int
+      get_array_length(const std::vector<types::global_dof_index> &array)
+      {
+        return array.size();
+      }
+
+      static unsigned int
+      get_array_length(const int &)
+      {
+        return 0;
+      }
+
+      static types::global_dof_index &
+      get_array_entry(std::vector<types::global_dof_index> &array,
+                      const unsigned int                    index)
+      {
+        AssertIndexRange(index, array.size());
+        return array[index];
+      }
+
+      static types::global_dof_index
+      get_array_entry(const std::vector<types::global_dof_index> &array,
+                      const unsigned int                          index)
+      {
+        AssertIndexRange(index, array.size());
+        return array[index];
+      }
+
+      static unsigned int
+      get_array_entry(const int &, const unsigned int)
+      {
+        return 0;
+      }
+
+
+
       /**
        * Process all dofs of an object.
        *
@@ -501,7 +543,7 @@ namespace internal
       template <int dim,
                 int spacedim,
                 int structdim,
-                typename GlobalIndexType,
+                typename DoFIndicesType,
                 typename DoFProcessor,
                 typename DoFMapping>
       static void
@@ -511,7 +553,7 @@ namespace internal
                      const unsigned int                            fe_index,
                      const DoFMapping &                            mapping,
                      const std::integral_constant<int, structdim> &dd,
-                     GlobalIndexType &   local_indices,
+                     DoFIndicesType &    local_indices,
                      unsigned int &      index,
                      const DoFProcessor &process)
       {
@@ -528,7 +570,7 @@ namespace internal
                     [structdim < dim ? 0 : obj_level][structdim]
                     [range.first +
                      ((structdim == 0 || structdim == dim) ? i : mapping(i))],
-                  local_indices[index]);
+                  get_array_entry(local_indices, index));
       }
 
 
@@ -597,39 +639,22 @@ namespace internal
 
 
       template <int dim, int spacedim>
-      static types::global_dof_index
-      mg_vertex_dof_index(const DoFHandler<dim, spacedim> &dof_handler,
-                          const int                        level,
-                          const unsigned int               vertex_index,
-                          const unsigned int               i)
+      static types::global_dof_index &
+      mg_vertex_dof_index(DoFHandler<dim, spacedim> &dof_handler,
+                          const int                  level,
+                          const unsigned int         vertex_index,
+                          const unsigned int         i)
       {
         Assert(dof_handler.hp_capability_enabled == false,
                ExcMessage(
                  "DoFHandler in hp-mode does not implement multilevel DoFs."));
 
-        return dof_handler.mg_vertex_dofs[vertex_index].get_index(
+        return dof_handler.mg_vertex_dofs[vertex_index].access_index(
           level, i, dof_handler.get_fe().n_dofs_per_vertex());
       }
 
 
 
-      template <int dim, int spacedim>
-      static void
-      set_mg_vertex_dof_index(DoFHandler<dim, spacedim> &dof_handler,
-                              const int                  level,
-                              const unsigned int         vertex_index,
-                              const unsigned int         i,
-                              types::global_dof_index    index)
-      {
-        Assert(dof_handler.hp_capability_enabled == false,
-               ExcMessage(
-                 "DoFHandler in hp-mode does not implement multilevel DoFs."));
-
-        return dof_handler.mg_vertex_dofs[vertex_index].set_index(
-          level, i, dof_handler.get_fe().n_dofs_per_vertex(), index);
-      }
-
-
       /**
        * Return the number of different finite elements that are active on a
        * given object such as a vertex, line, or cell.
@@ -982,23 +1007,26 @@ namespace internal
       static void
       process_dof_indices(
         const dealii::DoFAccessor<structdim, dim, spacedim, level_dof_access>
-          &                 accessor,
-        DoFIndicesType &    dof_indices,
-        const unsigned int  fe_index,
-        const DoFOperation &dof_operation,
-        const DoFProcessor &dof_processor,
-        const bool          count_level_dofs)
+          &                   accessor,
+        const DoFIndicesType &const_dof_indices,
+        const unsigned int    fe_index,
+        const DoFOperation &  dof_operation,
+        const DoFProcessor &  dof_processor,
+        const bool            count_level_dofs)
       {
         // we cannot rely on the template parameter level_dof_access here, since
         // the function get_mg_dof_indices()/set_mg_dof_indices() can be called
         // even if level_dof_access==false.
         (void)count_level_dofs;
 
-        AssertIndexRange(n_dof_indices(accessor, fe_index, count_level_dofs),
-                         dof_indices.size() + 1);
-
         const auto &fe = accessor.get_fe(fe_index);
 
+        // we want to pass in rvalue 'int' types as `DoFIndicesType`, but we
+        // need non-const references for std::vector<> types, so get in a
+        // const reference here and immediately cast it away
+        DoFIndicesType &dof_indices =
+          const_cast<DoFIndicesType &>(const_dof_indices);
+
         unsigned int index = 0;
 
         // 1) VERTEX dofs
@@ -1070,8 +1098,8 @@ namespace internal
         // get_dof_indices() also for artificial faces to determine if a face
         // is artificial.
         types::global_dof_index invalid_index = numbers::invalid_dof_index;
-        for (; index < dof_indices.size(); ++index)
-          dof_processor(invalid_index, dof_indices[index]);
+        for (; index < get_array_length(dof_indices); ++index)
+          dof_processor(invalid_index, get_array_entry(dof_indices, index));
       }
 
 
@@ -1086,14 +1114,14 @@ namespace internal
         /**
          * Return vertex DoF indices.
          */
-        template <int structdim, typename ArrayType, typename DoFProcessor>
+        template <int structdim, typename DoFIndicesType, typename DoFProcessor>
         DEAL_II_ALWAYS_INLINE void
         process_vertex_dofs(
           const dealii::DoFAccessor<structdim, dim, spacedim, level_dof_access>
             &                 accessor,
           const unsigned int  vertex,
           unsigned int &      index,
-          ArrayType &         index_value,
+          DoFIndicesType &    dof_indices,
           const unsigned int  fe_index_,
           const DoFProcessor &dof_processor) const
         {
@@ -1111,7 +1139,7 @@ namespace internal
               return d;
             },
             std::integral_constant<int, 0>(),
-            index_value,
+            dof_indices,
             index,
             dof_processor);
         }
@@ -1121,7 +1149,7 @@ namespace internal
          */
         template <int structdim,
                   typename DoFMapping,
-                  typename ArrayType,
+                  typename DoFIndicesType,
                   typename DoFProcessor>
         DEAL_II_ALWAYS_INLINE void
         process_dofs(
@@ -1129,7 +1157,7 @@ namespace internal
             &                 accessor,
           const DoFMapping &  mapping,
           unsigned int &      index,
-          ArrayType &         index_value,
+          DoFIndicesType &    dof_indices,
           const unsigned int  fe_index_,
           const DoFProcessor &dof_processor) const
         {
@@ -1143,7 +1171,7 @@ namespace internal
                          fe_index,
                          mapping,
                          std::integral_constant<int, structdim>(),
-                         index_value,
+                         dof_indices,
                          index,
                          dof_processor);
         }
@@ -1151,16 +1179,16 @@ namespace internal
         /**
          * Fallback for DoFInvalidAccessor.
          */
-        template <int structdim_,
+        template <int structdim,
                   typename DoFMapping,
-                  typename ArrayType,
+                  typename DoFIndicesType,
                   typename DoFProcessor>
         DEAL_II_ALWAYS_INLINE void
         process_dofs(
-          const dealii::DoFInvalidAccessor<structdim_, dim, spacedim> &,
+          const dealii::DoFInvalidAccessor<structdim, dim, spacedim> &,
           const DoFMapping &,
           unsigned int &,
-          ArrayType &,
+          DoFIndicesType &,
           const unsigned int,
           const DoFProcessor &) const
         {
@@ -1174,155 +1202,92 @@ namespace internal
        * An internal struct encapsulating the task of getting level (vertex)
        * DoF indices.
        */
-      template <int dim, int spacedim, bool level_dof_access, int structdim>
-      struct MGDoFIndexGetter
+      template <int dim, int spacedim>
+      struct MGDoFIndexProcessor
       {
         /**
          * Constructor.
          */
-        MGDoFIndexGetter(const FiniteElement<dim, spacedim> &fe,
-                         const unsigned int                  level)
-          : fe(fe)
-          , level(level)
+        MGDoFIndexProcessor(const unsigned int level)
+          : level(level)
         {}
 
         /**
          * Return vertex DoF indices.
          */
-        template <typename DoFProcessor>
+        template <int structdim,
+                  typename DoFIndicesType,
+                  typename DoFProcessor,
+                  bool level_dof_access>
         DEAL_II_ALWAYS_INLINE void
         process_vertex_dofs(
           const dealii::DoFAccessor<structdim, dim, spacedim, level_dof_access>
-            &                                   accessor,
-          const unsigned int                    vertex,
-          unsigned int &                        index,
-          std::vector<types::global_dof_index> &index_value,
-          const unsigned int                    fe_index,
-          const DoFProcessor &) const
-        {
-          (void)fe_index;
-
-          for (unsigned int d = 0; d < fe.template n_dofs_per_object<0>();
-               ++d, ++index)
-            index_value[index] = accessor.mg_vertex_dof_index(level, vertex, d);
-        }
-
-        /**
-         * Return DoF indices for lines, quads, and inner degrees of freedom.
-         */
-        template <int structdim_, typename DoFMapping, typename DoFProcessor>
-        DEAL_II_ALWAYS_INLINE void
-        process_dofs(
-          const dealii::DoFAccessor<structdim_, dim, spacedim, level_dof_access>
-            &                                   accessor,
-          const DoFMapping &                    mapping,
-          unsigned int &                        index,
-          std::vector<types::global_dof_index> &index_value,
-          const unsigned int                    fe_index,
-          const DoFProcessor &) const
+            &                 accessor,
+          const unsigned int  vertex,
+          unsigned int &      index,
+          DoFIndicesType &    dof_indices,
+          const unsigned int  fe_index,
+          const DoFProcessor &dof_processor) const
         {
-          (void)fe_index;
-
           for (unsigned int d = 0;
-               d < fe.template n_dofs_per_object<structdim_>(0);
+               d < accessor.get_fe(fe_index).template n_dofs_per_object<0>();
                ++d, ++index)
-            index_value[index] = accessor.mg_dof_index(level, mapping(d));
+            dof_processor(
+              accessor.dof_handler
+                ->mg_vertex_dofs[accessor.vertex_index(vertex)]
+                .access_index(level,
+                              d,
+                              accessor.get_fe(fe_index).n_dofs_per_vertex()),
+              get_array_entry(dof_indices, index));
         }
 
         /**
-         * Fallback for DoFInvalidAccessor.
+         * Return DoF indices for lines, quads, and inner degrees of freedom.
          */
-        template <int structdim_, typename DoFMapping, typename DoFProcessor>
+        template <int structdim,
+                  typename DoFMapping,
+                  typename DoFIndicesType,
+                  typename DoFProcessor,
+                  bool level_dof_access>
         DEAL_II_ALWAYS_INLINE void
         process_dofs(
-          const dealii::DoFInvalidAccessor<structdim_, dim, spacedim> &,
-          const DoFMapping &,
-          unsigned int &,
-          std::vector<types::global_dof_index> &,
-          const unsigned int,
-          const DoFProcessor &) const
-        {
-          Assert(false, ExcInternalError());
-        }
-
-      private:
-        const FiniteElement<dim, spacedim> &fe;
-        const unsigned int                  level;
-      };
-
-
-
-      /**
-       * An internal struct encapsulating the task of setting level (vertex)
-       * DoF indices.
-       */
-      template <int dim, int spacedim, bool level_dof_access, int structdim>
-      struct MGDoFIndexSetter
-      {
-        /**
-         * Constructor.
-         */
-        MGDoFIndexSetter(const FiniteElement<dim, spacedim> &fe,
-                         const unsigned int                  level)
-          : fe(fe)
-          , level(level)
-        {}
-
-        /**
-         * Set vertex DoF indices.
-         */
-        template <typename DoFProcessor>
-        DEAL_II_ALWAYS_INLINE void
-        process_vertex_dofs(
           const dealii::DoFAccessor<structdim, dim, spacedim, level_dof_access>
-            &                                         accessor,
-          const unsigned int                          vertex,
-          unsigned int &                              index,
-          const std::vector<types::global_dof_index> &index_value,
-          const unsigned int                          fe_index,
-          const DoFProcessor &) const
-        {
-          (void)fe_index;
-          for (unsigned int d = 0; d < fe.template n_dofs_per_object<0>();
-               ++d, ++index)
-            accessor.set_mg_vertex_dof_index(level,
-                                             vertex,
-                                             d,
-                                             index_value[index]);
-        }
-
-        /**
-         * Set DoF indices for lines, quads, and inner degrees of freedom.
-         */
-        template <int structdim_, typename DoFMapping, typename DoFProcessor>
-        DEAL_II_ALWAYS_INLINE void
-        process_dofs(
-          const dealii::DoFAccessor<structdim_, dim, spacedim, level_dof_access>
-            &                                         accessor,
-          const DoFMapping &                          mapping,
-          unsigned int &                              index,
-          const std::vector<types::global_dof_index> &index_value,
-          const unsigned int                          fe_index,
-          const DoFProcessor &) const
+            &                 accessor,
+          const DoFMapping &  mapping,
+          unsigned int &      index,
+          DoFIndicesType &    dof_indices,
+          const unsigned int  fe_index,
+          const DoFProcessor &dof_processor) const
         {
-          (void)fe_index;
-
           for (unsigned int d = 0;
-               d < fe.template n_dofs_per_object<structdim_>(0);
+               d <
+               accessor.get_fe(fe_index).template n_dofs_per_object<structdim>(
+                 0);
                ++d, ++index)
-            accessor.set_mg_dof_index(level, mapping(d), index_value[index]);
+            dof_processor(
+              get_mg_dof_index(*accessor.dof_handler,
+                               accessor.dof_handler->mg_levels[level],
+                               accessor.dof_handler->mg_faces,
+                               accessor.index(),
+                               fe_index,
+                               (structdim == dim ? d : mapping(d)),
+                               std::integral_constant<int, structdim>()),
+              get_array_entry(dof_indices, index));
         }
 
         /**
          * Fallback for DoFInvalidAccessor.
          */
-        template <int structdim_, typename DoFMapping, typename DoFProcessor>
+        template <int structdim,
+                  typename DoFMapping,
+                  typename DoFIndicesType,
+                  typename DoFProcessor>
         DEAL_II_ALWAYS_INLINE void
         process_dofs(
-          const dealii::DoFInvalidAccessor<structdim_, dim, spacedim> &,
+          const dealii::DoFInvalidAccessor<structdim, dim, spacedim> &,
           const DoFMapping &,
           unsigned int &,
-          const std::vector<types::global_dof_index> &,
+          DoFIndicesType &,
           const unsigned int,
           const DoFProcessor &) const
         {
@@ -1330,8 +1295,7 @@ namespace internal
         }
 
       private:
-        const FiniteElement<dim, spacedim> &fe;
-        const unsigned int                  level;
+        const unsigned int level;
       };
 
 
@@ -1349,7 +1313,7 @@ namespace internal
           dof_indices,
           fe_index,
           DoFIndexProcessor<dim, spacedim, level_dof_access>(),
-          [](auto &stored_index, auto &dof_index) { dof_index = stored_index; },
+          [](auto stored_index, auto &dof_index) { dof_index = stored_index; },
           false);
       }
 
@@ -1379,7 +1343,7 @@ namespace internal
           dof_indices,
           fe_index,
           DoFIndexProcessor<dim, spacedim, level_dof_access>(),
-          [](auto &stored_index, auto &dof_index) { stored_index = dof_index; },
+          [](auto &stored_index, auto dof_index) { stored_index = dof_index; },
           false);
       }
 
@@ -1398,9 +1362,8 @@ namespace internal
           accessor,
           dof_indices,
           fe_index,
-          MGDoFIndexGetter<dim, spacedim, level_dof_access, structdim>(
-            accessor.get_fe(fe_index), level),
-          [](auto &, auto &) {},
+          MGDoFIndexProcessor<dim, spacedim>(level),
+          [](auto stored_index, auto &dof_index) { dof_index = stored_index; },
           true);
       }
 
@@ -1420,20 +1383,89 @@ namespace internal
         // function DoFCellAccessor::set_mg_dof_indices(), which is called by
         // internal::DoFHandlerImplementation::Policy::Implementation::distribute_mg_dofs().
         // In the case of new use cases, this assert can be removed.
-        Assert(
-          dim == structdim,
-          ExcMessage(
-            "This function is intended to be used for DoFCellAccessor, i.e., dimension == structdim."));
+        Assert(dim == structdim,
+               ExcMessage("This function is intended to be used for "
+                          "DoFCellAccessor, i.e., dimension == structdim."));
 
         process_dof_indices(
           accessor,
           dof_indices,
           fe_index,
-          MGDoFIndexSetter<dim, spacedim, level_dof_access, structdim>(
-            accessor.get_fe(fe_index), level),
-          [](auto &, auto &) {},
+          MGDoFIndexProcessor<dim, spacedim>(level),
+          [](auto &stored_index, auto dof_index) { stored_index = dof_index; },
           true);
       }
+
+
+
+      template <int dim, int spacedim>
+      static types::global_dof_index &
+      get_mg_dof_index(
+        const DoFHandler<dim, spacedim> &dof_handler,
+        const std::unique_ptr<internal::DoFHandlerImplementation::DoFLevel<dim>>
+          &mg_level,
+        const std::unique_ptr<internal::DoFHandlerImplementation::DoFFaces<dim>>
+          &,
+        const unsigned int obj_index,
+        const unsigned int fe_index,
+        const unsigned int local_index,
+        const std::integral_constant<int, dim>)
+      {
+        Assert(dof_handler.hp_capability_enabled == false,
+               (typename DoFHandler<dim, spacedim>::ExcNotImplementedWithHP()));
+
+        return mg_level->dof_object.access_dof_index(
+          static_cast<const DoFHandler<dim, spacedim> &>(dof_handler),
+          obj_index,
+          fe_index,
+          local_index);
+      }
+
+
+
+      template <int dim, int spacedim, std::enable_if_t<(dim > 1), int> = 0>
+      static types::global_dof_index &
+      get_mg_dof_index(
+        const DoFHandler<dim, spacedim> &dof_handler,
+        const std::unique_ptr<internal::DoFHandlerImplementation::DoFLevel<dim>>
+          &,
+        const std::unique_ptr<internal::DoFHandlerImplementation::DoFFaces<dim>>
+          &                mg_faces,
+        const unsigned int obj_index,
+        const unsigned int fe_index,
+        const unsigned int local_index,
+        const std::integral_constant<int, 1>)
+      {
+        return mg_faces->lines.access_dof_index(
+          static_cast<const DoFHandler<dim, spacedim> &>(dof_handler),
+          obj_index,
+          fe_index,
+          local_index);
+      }
+
+
+
+      template <int spacedim>
+      static types::global_dof_index &
+      get_mg_dof_index(
+        const DoFHandler<3, spacedim> &dof_handler,
+        const std::unique_ptr<internal::DoFHandlerImplementation::DoFLevel<3>>
+          &,
+        const std::unique_ptr<internal::DoFHandlerImplementation::DoFFaces<3>>
+          &                mg_faces,
+        const unsigned int obj_index,
+        const unsigned int fe_index,
+        const unsigned int local_index,
+        const std::integral_constant<int, 2>)
+      {
+        Assert(dof_handler.hp_capability_enabled == false,
+               (typename DoFHandler<3, spacedim>::ExcNotImplementedWithHP()));
+        return mg_faces->quads.access_dof_index(
+          static_cast<const DoFHandler<3, spacedim> &>(dof_handler),
+          obj_index,
+          fe_index,
+          local_index);
+      }
     };
   } // namespace DoFAccessorImplementation
 } // namespace internal
@@ -1467,10 +1499,14 @@ DoFAccessor<structdim, dim, spacedim, level_dof_access>::mg_dof_index(
   const int          level,
   const unsigned int i) const
 {
-  return this->dof_handler->template get_dof_index<structdim>(level,
-                                                              this->index(),
-                                                              0,
-                                                              i);
+  return internal::DoFAccessorImplementation::Implementation::get_mg_dof_index(
+    *this->dof_handler,
+    this->dof_handler->mg_levels[level],
+    this->dof_handler->mg_faces,
+    this->index(),
+    0,
+    i,
+    std::integral_constant<int, structdim>());
 }
 
 
@@ -1498,6 +1534,25 @@ DoFAccessor<structdim, dim, spacedim, level_dof_access>::set_dof_index(
 
 
 
+template <int structdim, int dim, int spacedim, bool level_dof_access>
+inline void
+DoFAccessor<structdim, dim, spacedim, level_dof_access>::set_mg_dof_index(
+  const int                     level,
+  const unsigned int            i,
+  const types::global_dof_index index) const
+{
+  internal::DoFAccessorImplementation::Implementation::get_mg_dof_index(
+    *this->dof_handler,
+    this->dof_handler->mg_levels[level],
+    this->dof_handler->mg_faces,
+    this->index(),
+    0,
+    i,
+    std::integral_constant<int, structdim>()) = index;
+}
+
+
+
 template <int structdim, int dim, int spacedim, bool level_dof_access>
 inline unsigned int
 DoFAccessor<structdim, dim, spacedim, level_dof_access>::n_active_fe_indices()
@@ -1615,7 +1670,7 @@ DoFAccessor<structdim, dim, spacedim, level_dof_access>::mg_vertex_dof_index(
            "DoFHandler in hp-mode does not implement multilevel DoFs."));
 
   return this->dof_handler->mg_vertex_dofs[this->vertex_index(vertex)]
-    .get_index(level, i, this->dof_handler->get_fe().n_dofs_per_vertex());
+    .access_index(level, i, this->dof_handler->get_fe().n_dofs_per_vertex());
 }
 
 
@@ -1641,20 +1696,8 @@ DoFAccessor<structdim, dim, spacedim, level_dof_access>::
          ExcMessage(
            "DoFHandler in hp-mode does not implement multilevel DoFs."));
 
-  this->dof_handler->mg_vertex_dofs[this->vertex_index(vertex)].set_index(
-    level, i, this->dof_handler->get_fe().n_dofs_per_vertex(), index);
-}
-
-
-template <int structdim, int dim, int spacedim, bool level_dof_access>
-inline void
-DoFAccessor<structdim, dim, spacedim, level_dof_access>::set_mg_dof_index(
-  const int                     level,
-  const unsigned int            i,
-  const types::global_dof_index index) const
-{
-  this->dof_handler->template set_dof_index<structdim>(
-    level, this->index(), 0, i, index);
+  this->dof_handler->mg_vertex_dofs[this->vertex_index(vertex)].access_index(
+    level, i, this->dof_handler->get_fe().n_dofs_per_vertex()) = index;
 }
 
 
index 12f74742eda6d37b8d27c50bc22e081468561541..f850febf210e12192479e147b52367a789a5cd8b 100644 (file)
@@ -1391,20 +1391,10 @@ private:
      * Return the index of the <code>dof_number</code>th degree of freedom for
      * the given level stored for the current vertex.
      */
-    types::global_dof_index
-    get_index(const unsigned int level,
-              const unsigned int dof_number,
-              const unsigned int dofs_per_vertex) const;
-
-    /**
-     * Set the index of the <code>dof_number</code>th degree of freedom for
-     * the given level stored for the current vertex to <code>index</code>.
-     */
-    void
-    set_index(const unsigned int            level,
-              const unsigned int            dof_number,
-              const unsigned int            dofs_per_vertex,
-              const types::global_dof_index index);
+    types::global_dof_index &
+    access_index(const unsigned int level,
+                 const unsigned int dof_number,
+                 const unsigned int dofs_per_vertex);
 
   private:
     /**
@@ -1630,27 +1620,6 @@ private:
   void
   clear_mg_space();
 
-  /**
-   * Return dof index of specified object.
-   */
-  template <int structdim>
-  types::global_dof_index
-  get_dof_index(const unsigned int obj_level,
-                const unsigned int obj_index,
-                const unsigned int fe_index,
-                const unsigned int local_index) const;
-
-  /**
-   * Return dof index of specified object.
-   */
-  template <int structdim>
-  void
-  set_dof_index(const unsigned int            obj_level,
-                const unsigned int            obj_index,
-                const unsigned int            fe_index,
-                const unsigned int            local_index,
-                const types::global_dof_index global_index) const;
-
   /**
    * Set up DoFHandler policy.
    */
@@ -2156,11 +2125,11 @@ DoFHandler<dim, spacedim>::load(Archive &ar, const unsigned int)
 
 
 template <int dim, int spacedim>
-inline types::global_dof_index
-DoFHandler<dim, spacedim>::MGVertexDoFs::get_index(
+inline types::global_dof_index &
+DoFHandler<dim, spacedim>::MGVertexDoFs::access_index(
   const unsigned int level,
   const unsigned int dof_number,
-  const unsigned int dofs_per_vertex) const
+  const unsigned int dofs_per_vertex)
 {
   Assert((level >= coarsest_level) && (level <= finest_level),
          ExcInvalidLevel(level));
@@ -2169,21 +2138,6 @@ DoFHandler<dim, spacedim>::MGVertexDoFs::get_index(
 
 
 
-template <int dim, int spacedim>
-inline void
-DoFHandler<dim, spacedim>::MGVertexDoFs::set_index(
-  const unsigned int            level,
-  const unsigned int            dof_number,
-  const unsigned int            dofs_per_vertex,
-  const types::global_dof_index index)
-{
-  Assert((level >= coarsest_level) && (level <= finest_level),
-         ExcInvalidLevel(level));
-  indices[dofs_per_vertex * (level - coarsest_level) + dof_number] = index;
-}
-
-
-
 extern template class DoFHandler<1, 1>;
 extern template class DoFHandler<1, 2>;
 extern template class DoFHandler<1, 3>;
index b563649d68bd610d98e357b14a45b1f034311946..151ba7551b80ccfa39fe2414a5d9462effb080ba 100644 (file)
@@ -74,28 +74,6 @@ namespace internal
        */
       std::vector<types::global_dof_index> dofs;
 
-    public:
-      /**
-       * Set the global index of the @p local_index-th degree of freedom
-       * located on the object with number @p obj_index to the value given by
-       * the last argument. The @p dof_handler argument is used to access the
-       * finite element that is to be used to compute the location where this
-       * data is stored.
-       *
-       * The third argument, @p fe_index, must equal zero. It is otherwise
-       * unused, but we retain the argument so that we can use the same
-       * interface for non-hp- and hp-finite element methods, in effect making
-       * it possible to share the DoFAccessor class hierarchy between hp- and
-       * non-hp-classes.
-       */
-      template <int dh_dim, int spacedim>
-      void
-      set_dof_index(const dealii::DoFHandler<dh_dim, spacedim> &dof_handler,
-                    const unsigned int                          obj_index,
-                    const unsigned int                          fe_index,
-                    const unsigned int                          local_index,
-                    const types::global_dof_index               global_index);
-
       /**
        * Return the global index of the @p local_index-th degree of freedom
        * located on the object with number @p obj_index. The @p dof_handler
@@ -109,11 +87,11 @@ namespace internal
        * non-hp-classes.
        */
       template <int dh_dim, int spacedim>
-      types::global_dof_index
-      get_dof_index(const dealii::DoFHandler<dh_dim, spacedim> &dof_handler,
-                    const unsigned int                          obj_index,
-                    const unsigned int                          fe_index,
-                    const unsigned int local_index) const;
+      types::global_dof_index &
+      access_dof_index(const dealii::DoFHandler<dh_dim, spacedim> &dof_handler,
+                       const unsigned int                          obj_index,
+                       const unsigned int                          fe_index,
+                       const unsigned int                          local_index);
 
       /**
        * Return the value 1. The meaning of this function becomes clear by
@@ -196,12 +174,12 @@ namespace internal
 
     template <int dim>
     template <int dh_dim, int spacedim>
-    inline types::global_dof_index
-    DoFObjects<dim>::get_dof_index(
+    inline types::global_dof_index &
+    DoFObjects<dim>::access_dof_index(
       const dealii::DoFHandler<dh_dim, spacedim> &dof_handler,
       const unsigned int                          obj_index,
       const unsigned int                          fe_index,
-      const unsigned int                          local_index) const
+      const unsigned int                          local_index)
     {
       (void)fe_index;
       Assert(
index 33f7ccd18d00c9b11fad6641dcc8166b468280cd..c7aaa49023d02e18df4e2fdb4515136a9ac27b66 100644 (file)
@@ -673,281 +673,6 @@ namespace internal
               dof_handler.mg_vertex_dofs[vertex].init(1, 0, 0);
             }
       }
-
-      template <int spacedim>
-      static types::global_dof_index
-      get_dof_index(
-        const DoFHandler<1, spacedim> &dof_handler,
-        const std::unique_ptr<internal::DoFHandlerImplementation::DoFLevel<1>>
-          &mg_level,
-        const std::unique_ptr<internal::DoFHandlerImplementation::DoFFaces<1>>
-          &,
-        const unsigned int obj_index,
-        const unsigned int fe_index,
-        const unsigned int local_index,
-        const std::integral_constant<int, 1>)
-      {
-        Assert(dof_handler.hp_capability_enabled == false,
-               (typename DoFHandler<1, spacedim>::ExcNotImplementedWithHP()));
-
-        return mg_level->dof_object.get_dof_index(
-          static_cast<const DoFHandler<1, spacedim> &>(dof_handler),
-          obj_index,
-          fe_index,
-          local_index);
-      }
-
-      template <int spacedim>
-      static types::global_dof_index
-      get_dof_index(
-        const DoFHandler<2, spacedim> &dof_handler,
-        const std::unique_ptr<internal::DoFHandlerImplementation::DoFLevel<2>>
-          &,
-        const std::unique_ptr<internal::DoFHandlerImplementation::DoFFaces<2>>
-          &                mg_faces,
-        const unsigned int obj_index,
-        const unsigned int fe_index,
-        const unsigned int local_index,
-        const std::integral_constant<int, 1>)
-      {
-        return mg_faces->lines.get_dof_index(
-          static_cast<const DoFHandler<2, spacedim> &>(dof_handler),
-          obj_index,
-          fe_index,
-          local_index);
-      }
-
-      template <int spacedim>
-      static types::global_dof_index
-      get_dof_index(
-        const DoFHandler<2, spacedim> &dof_handler,
-        const std::unique_ptr<internal::DoFHandlerImplementation::DoFLevel<2>>
-          &mg_level,
-        const std::unique_ptr<internal::DoFHandlerImplementation::DoFFaces<2>>
-          &,
-        const unsigned int obj_index,
-        const unsigned int fe_index,
-        const unsigned int local_index,
-        const std::integral_constant<int, 2>)
-      {
-        Assert(dof_handler.hp_capability_enabled == false,
-               (typename DoFHandler<2, spacedim>::ExcNotImplementedWithHP()));
-        return mg_level->dof_object.get_dof_index(
-          static_cast<const DoFHandler<2, spacedim> &>(dof_handler),
-          obj_index,
-          fe_index,
-          local_index);
-      }
-
-      template <int spacedim>
-      static types::global_dof_index
-      get_dof_index(
-        const DoFHandler<3, spacedim> &dof_handler,
-        const std::unique_ptr<internal::DoFHandlerImplementation::DoFLevel<3>>
-          &,
-        const std::unique_ptr<internal::DoFHandlerImplementation::DoFFaces<3>>
-          &                mg_faces,
-        const unsigned int obj_index,
-        const unsigned int fe_index,
-        const unsigned int local_index,
-        const std::integral_constant<int, 1>)
-      {
-        Assert(dof_handler.hp_capability_enabled == false,
-               (typename DoFHandler<3, spacedim>::ExcNotImplementedWithHP()));
-        return mg_faces->lines.get_dof_index(
-          static_cast<const DoFHandler<3, spacedim> &>(dof_handler),
-          obj_index,
-          fe_index,
-          local_index);
-      }
-
-      template <int spacedim>
-      static types::global_dof_index
-      get_dof_index(
-        const DoFHandler<3, spacedim> &dof_handler,
-        const std::unique_ptr<internal::DoFHandlerImplementation::DoFLevel<3>>
-          &,
-        const std::unique_ptr<internal::DoFHandlerImplementation::DoFFaces<3>>
-          &                mg_faces,
-        const unsigned int obj_index,
-        const unsigned int fe_index,
-        const unsigned int local_index,
-        const std::integral_constant<int, 2>)
-      {
-        Assert(dof_handler.hp_capability_enabled == false,
-               (typename DoFHandler<3, spacedim>::ExcNotImplementedWithHP()));
-        return mg_faces->quads.get_dof_index(
-          static_cast<const DoFHandler<3, spacedim> &>(dof_handler),
-          obj_index,
-          fe_index,
-          local_index);
-      }
-
-      template <int spacedim>
-      static types::global_dof_index
-      get_dof_index(
-        const DoFHandler<3, spacedim> &dof_handler,
-        const std::unique_ptr<internal::DoFHandlerImplementation::DoFLevel<3>>
-          &mg_level,
-        const std::unique_ptr<internal::DoFHandlerImplementation::DoFFaces<3>>
-          &,
-        const unsigned int obj_index,
-        const unsigned int fe_index,
-        const unsigned int local_index,
-        const std::integral_constant<int, 3>)
-      {
-        Assert(dof_handler.hp_capability_enabled == false,
-               (typename DoFHandler<3, spacedim>::ExcNotImplementedWithHP()));
-        return mg_level->dof_object.get_dof_index(
-          static_cast<const DoFHandler<3, spacedim> &>(dof_handler),
-          obj_index,
-          fe_index,
-          local_index);
-      }
-
-      template <int spacedim>
-      static void
-      set_dof_index(
-        const DoFHandler<1, spacedim> &dof_handler,
-        const std::unique_ptr<internal::DoFHandlerImplementation::DoFLevel<1>>
-          &mg_level,
-        const std::unique_ptr<internal::DoFHandlerImplementation::DoFFaces<1>>
-          &,
-        const unsigned int            obj_index,
-        const unsigned int            fe_index,
-        const unsigned int            local_index,
-        const types::global_dof_index global_index,
-        const std::integral_constant<int, 1>)
-      {
-        Assert(dof_handler.hp_capability_enabled == false,
-               (typename DoFHandler<1, spacedim>::ExcNotImplementedWithHP()));
-        mg_level->dof_object.set_dof_index(
-          static_cast<const DoFHandler<1, spacedim> &>(dof_handler),
-          obj_index,
-          fe_index,
-          local_index,
-          global_index);
-      }
-
-      template <int spacedim>
-      static void
-      set_dof_index(
-        const DoFHandler<2, spacedim> &dof_handler,
-        const std::unique_ptr<internal::DoFHandlerImplementation::DoFLevel<2>>
-          &,
-        const std::unique_ptr<internal::DoFHandlerImplementation::DoFFaces<2>>
-          &                           mg_faces,
-        const unsigned int            obj_index,
-        const unsigned int            fe_index,
-        const unsigned int            local_index,
-        const types::global_dof_index global_index,
-        const std::integral_constant<int, 1>)
-      {
-        Assert(dof_handler.hp_capability_enabled == false,
-               (typename DoFHandler<2, spacedim>::ExcNotImplementedWithHP()));
-        mg_faces->lines.set_dof_index(
-          static_cast<const DoFHandler<2, spacedim> &>(dof_handler),
-          obj_index,
-          fe_index,
-          local_index,
-          global_index);
-      }
-
-      template <int spacedim>
-      static void
-      set_dof_index(
-        const DoFHandler<2, spacedim> &dof_handler,
-        const std::unique_ptr<internal::DoFHandlerImplementation::DoFLevel<2>>
-          &mg_level,
-        const std::unique_ptr<internal::DoFHandlerImplementation::DoFFaces<2>>
-          &,
-        const unsigned int            obj_index,
-        const unsigned int            fe_index,
-        const unsigned int            local_index,
-        const types::global_dof_index global_index,
-        const std::integral_constant<int, 2>)
-      {
-        Assert(dof_handler.hp_capability_enabled == false,
-               (typename DoFHandler<2, spacedim>::ExcNotImplementedWithHP()));
-        mg_level->dof_object.set_dof_index(
-          static_cast<const DoFHandler<2, spacedim> &>(dof_handler),
-          obj_index,
-          fe_index,
-          local_index,
-          global_index);
-      }
-
-      template <int spacedim>
-      static void
-      set_dof_index(
-        const DoFHandler<3, spacedim> &dof_handler,
-        const std::unique_ptr<internal::DoFHandlerImplementation::DoFLevel<3>>
-          &,
-        const std::unique_ptr<internal::DoFHandlerImplementation::DoFFaces<3>>
-          &                           mg_faces,
-        const unsigned int            obj_index,
-        const unsigned int            fe_index,
-        const unsigned int            local_index,
-        const types::global_dof_index global_index,
-        const std::integral_constant<int, 1>)
-      {
-        Assert(dof_handler.hp_capability_enabled == false,
-               (typename DoFHandler<3, spacedim>::ExcNotImplementedWithHP()));
-        mg_faces->lines.set_dof_index(
-          static_cast<const DoFHandler<3, spacedim> &>(dof_handler),
-          obj_index,
-          fe_index,
-          local_index,
-          global_index);
-      }
-
-      template <int spacedim>
-      static void
-      set_dof_index(
-        const DoFHandler<3, spacedim> &dof_handler,
-        const std::unique_ptr<internal::DoFHandlerImplementation::DoFLevel<3>>
-          &,
-        const std::unique_ptr<internal::DoFHandlerImplementation::DoFFaces<3>>
-          &                           mg_faces,
-        const unsigned int            obj_index,
-        const unsigned int            fe_index,
-        const unsigned int            local_index,
-        const types::global_dof_index global_index,
-        const std::integral_constant<int, 2>)
-      {
-        Assert(dof_handler.hp_capability_enabled == false,
-               (typename DoFHandler<3, spacedim>::ExcNotImplementedWithHP()));
-        mg_faces->quads.set_dof_index(
-          static_cast<const DoFHandler<3, spacedim> &>(dof_handler),
-          obj_index,
-          fe_index,
-          local_index,
-          global_index);
-      }
-
-      template <int spacedim>
-      static void
-      set_dof_index(
-        const DoFHandler<3, spacedim> &dof_handler,
-        const std::unique_ptr<internal::DoFHandlerImplementation::DoFLevel<3>>
-          &mg_level,
-        const std::unique_ptr<internal::DoFHandlerImplementation::DoFFaces<3>>
-          &,
-        const unsigned int            obj_index,
-        const unsigned int            fe_index,
-        const unsigned int            local_index,
-        const types::global_dof_index global_index,
-        const std::integral_constant<int, 3>)
-      {
-        Assert(dof_handler.hp_capability_enabled == false,
-               (typename DoFHandler<3, spacedim>::ExcNotImplementedWithHP()));
-        mg_level->dof_object.set_dof_index(
-          static_cast<const DoFHandler<3, spacedim> &>(dof_handler),
-          obj_index,
-          fe_index,
-          local_index,
-          global_index);
-      }
     };
   } // namespace DoFHandlerImplementation
 
@@ -3069,65 +2794,6 @@ DoFHandler<dim, spacedim>::max_couplings_between_dofs() const
 
 
 
-template <int dim, int spacedim>
-template <int structdim>
-types::global_dof_index
-DoFHandler<dim, spacedim>::get_dof_index(const unsigned int obj_level,
-                                         const unsigned int obj_index,
-                                         const unsigned int fe_index,
-                                         const unsigned int local_index) const
-{
-  if (hp_capability_enabled)
-    {
-      Assert(false, ExcNotImplemented());
-      return numbers::invalid_dof_index;
-    }
-  else
-    {
-      return internal::DoFHandlerImplementation::Implementation::get_dof_index(
-        *this,
-        this->mg_levels[obj_level],
-        this->mg_faces,
-        obj_index,
-        fe_index,
-        local_index,
-        std::integral_constant<int, structdim>());
-    }
-}
-
-
-
-template <int dim, int spacedim>
-template <int structdim>
-void
-DoFHandler<dim, spacedim>::set_dof_index(
-  const unsigned int            obj_level,
-  const unsigned int            obj_index,
-  const unsigned int            fe_index,
-  const unsigned int            local_index,
-  const types::global_dof_index global_index) const
-{
-  if (hp_capability_enabled)
-    {
-      Assert(false, ExcNotImplemented());
-      return;
-    }
-  else
-    {
-      internal::DoFHandlerImplementation::Implementation::set_dof_index(
-        *this,
-        this->mg_levels[obj_level],
-        this->mg_faces,
-        obj_index,
-        fe_index,
-        local_index,
-        global_index,
-        std::integral_constant<int, structdim>());
-    }
-}
-
-
-
 template <int dim, int spacedim>
 void
 DoFHandler<dim, spacedim>::set_active_fe_indices(
index e78524dd3102fd5d2bdb6fb9704e3b1685cda525..59075cffe894c54eaa11a199570f243e2da6000f 100644 (file)
@@ -46,26 +46,3 @@ for (deal_II_dimension : DIMENSIONS; deal_II_space_dimension : DIMENSIONS)
     \}
 #endif
   }
-
-
-for (deal_II_dimension : DIMENSIONS; deal_II_space_dimension : DIMENSIONS;
-     deal_II_component : DIMENSIONS)
-  {
-#if deal_II_component <= deal_II_dimension && \
-  deal_II_dimension <= deal_II_space_dimension
-    template types::global_dof_index
-    DoFHandler<deal_II_dimension, deal_II_space_dimension>::get_dof_index<
-      deal_II_component>(const unsigned int obj_level,
-                         const unsigned int obj_index,
-                         const unsigned int fe_index,
-                         const unsigned int local_index) const;
-
-    template void
-    DoFHandler<deal_II_dimension, deal_II_space_dimension>::set_dof_index<
-      deal_II_component>(const unsigned int            obj_level,
-                         const unsigned int            obj_index,
-                         const unsigned int            fe_index,
-                         const unsigned int            local_index,
-                         const types::global_dof_index global_index) const;
-#endif
-  }
index d9de7cb4da902e1d89f204a35f18912bf07f8dfe..c8ab7b57b2ec2cac0279a97707a02ee44ae233b3 100644 (file)
@@ -1758,22 +1758,21 @@ namespace internal
           // distribute dofs on all cells excluding artificial ones
           types::global_dof_index next_free_dof = 0;
 
-          std::vector<types::global_dof_index> dof_indices;
-
           for (auto cell : dof_handler.active_cell_iterators())
             if (!cell->is_artificial() &&
                 ((subdomain_id == numbers::invalid_subdomain_id) ||
                  (cell->subdomain_id() == subdomain_id)))
               {
-                dof_indices.resize(cell->get_fe().n_dofs_per_cell());
-
+                // need to create a dummy int reference to feed the
+                // process_dof_indices function for the case we do want want
+                // to retrieve any index
                 DoFAccessorImplementation::Implementation::process_dof_indices(
                   *cell,
-                  dof_indices,
+                  0,
                   cell->active_fe_index(),
                   DoFAccessorImplementation::Implementation::
                     DoFIndexProcessor<dim, spacedim, false>(),
-                  [&next_free_dof](auto &stored_index, auto &) {
+                  [&next_free_dof](auto &stored_index, auto) {
                     if (stored_index == numbers::invalid_dof_index)
                       {
                         stored_index = next_free_dof;
@@ -1863,35 +1862,34 @@ namespace internal
 
           types::global_dof_index next_free_dof = 0;
 
-          std::vector<types::global_dof_index> dof_indices;
-
           for (auto cell : dof_handler.cell_iterators_on_level(level))
             if ((level_subdomain_id == numbers::invalid_subdomain_id) ||
                 (cell->level_subdomain_id() == level_subdomain_id))
               {
-                dof_indices.resize(cell->get_fe().n_dofs_per_cell());
-
-                cell->get_mg_dof_indices(dof_indices);
-
-                for (auto &dof_index : dof_indices)
-                  if (dof_index == numbers::invalid_dof_index)
-                    {
-                      dof_index = next_free_dof;
-
-                      Assert(
-                        next_free_dof !=
-                          std::numeric_limits<types::global_dof_index>::max(),
-                        ExcMessage(
-                          "You have reached the maximal number of degrees of "
-                          "freedom that can be stored in the chosen data type. "
-                          "In practice, this can only happen if you are using "
-                          "32-bit data types. You will have to re-compile "
-                          "deal.II with the `DEAL_II_WITH_64BIT_INDICES' "
-                          "flag set to `ON'."));
-                      ++next_free_dof;
-                    }
-
-                cell->set_mg_dof_indices(dof_indices);
+                DoFAccessorImplementation::Implementation::process_dof_indices(
+                  *cell,
+                  0,
+                  0,
+                  DoFAccessorImplementation::Implementation::
+                    MGDoFIndexProcessor<dim, spacedim>(level),
+                  [&next_free_dof](auto &stored_index, auto) {
+                    if (stored_index == numbers::invalid_dof_index)
+                      {
+                        stored_index = next_free_dof;
+                        Assert(
+                          next_free_dof !=
+                            std::numeric_limits<types::global_dof_index>::max(),
+                          ExcMessage(
+                            "You have reached the maximal number of degrees of "
+                            "freedom that can be stored in the chosen data "
+                            "type. In practice, this can only happen if you "
+                            "are using 32-bit data types. You will have to "
+                            "re-compile deal.II with the "
+                            "`DEAL_II_WITH_64BIT_INDICES' flag set to `ON'."));
+                        ++next_free_dof;
+                      }
+                  },
+                  true);
               }
 
           return next_free_dof;
@@ -2549,9 +2547,9 @@ namespace internal
                    ++d)
                 {
                   const dealii::types::global_dof_index idx =
-                    i->get_index(level,
-                                 d,
-                                 dof_handler.get_fe().n_dofs_per_vertex());
+                    i->access_index(level,
+                                    d,
+                                    dof_handler.get_fe().n_dofs_per_vertex());
 
                   if (idx != numbers::invalid_dof_index)
                     {
@@ -2559,13 +2557,12 @@ namespace internal
                                indices_we_care_about.is_element(idx) :
                                (idx < new_numbers.size()),
                              ExcInternalError());
-                      i->set_index(level,
-                                   d,
-                                   dof_handler.get_fe().n_dofs_per_vertex(),
-                                   (indices_we_care_about.size() == 0) ?
-                                     (new_numbers[idx]) :
-                                     (new_numbers[indices_we_care_about
-                                                    .index_within_set(idx)]));
+                      i->access_index(
+                        level, d, dof_handler.get_fe().n_dofs_per_vertex()) =
+                        (indices_we_care_about.size() == 0) ?
+                          new_numbers[idx] :
+                          new_numbers[indices_we_care_about.index_within_set(
+                            idx)];
                     }
                 }
         }
@@ -3634,21 +3631,27 @@ namespace internal
                      dealii::numbers::artificial_subdomain_id,
                    ExcInternalError());
 
-            std::vector<dealii::types::global_dof_index> dof_indices(
-              cell->get_fe().n_dofs_per_cell());
-            cell->get_mg_dof_indices(dof_indices);
-
             bool complete = true;
-            for (unsigned int i = 0; i < dof_indices.size(); ++i)
-              if (dofs[i] != numbers::invalid_dof_index)
-                {
-                  Assert((dof_indices[i] == (numbers::invalid_dof_index)) ||
-                           (dof_indices[i] == dofs[i]),
-                         ExcInternalError());
-                  dof_indices[i] = dofs[i];
-                }
-              else
-                complete = false;
+            DoFAccessorImplementation::Implementation::process_dof_indices(
+              *cell,
+              dofs,
+              0,
+              DoFAccessorImplementation::Implementation::
+                MGDoFIndexProcessor<dim, spacedim>(cell->level()),
+              [&complete](auto &stored_index, auto received_index) {
+                if (received_index != numbers::invalid_dof_index)
+                  {
+#  if !defined(__INTEL_COMPILER) || __INTEL_COMPILER >= 1900
+                    Assert((stored_index == (numbers::invalid_dof_index)) ||
+                             (stored_index == received_index),
+                           ExcInternalError());
+#  endif
+                    stored_index = received_index;
+                  }
+                else
+                  complete = false;
+              },
+              true);
 
             if (!complete)
               const_cast<
@@ -3658,10 +3661,6 @@ namespace internal
               const_cast<
                 typename DoFHandler<dim, spacedim>::level_cell_iterator &>(cell)
                 ->clear_user_flag();
-
-            const_cast<
-              typename DoFHandler<dim, spacedim>::level_cell_iterator &>(cell)
-              ->set_mg_dof_indices(dof_indices);
           };
 
           const auto filter = [](const auto &cell) {
@@ -3735,7 +3734,7 @@ namespace internal
               cell->active_fe_index(),
               DoFAccessorImplementation::Implementation::
                 DoFIndexProcessor<dim, spacedim, false>(),
-              [&complete](auto &stored_index, auto &received_index) {
+              [&complete](auto &stored_index, auto received_index) {
                 if (received_index != numbers::invalid_dof_index)
                   {
 #    if !defined(__INTEL_COMPILER) || __INTEL_COMPILER >= 1900
@@ -4272,26 +4271,25 @@ namespace internal
         // locally owned cell and a ghost cell. In any case, it is
         // sufficient to kill them only from the ghost side cell, so loop
         // only over ghost cells
-        {
-          std::vector<dealii::types::global_dof_index> local_dof_indices;
-
-          for (auto cell : dof_handler->active_cell_iterators())
-            if (cell->is_ghost())
-              {
-                local_dof_indices.resize(cell->get_fe().n_dofs_per_cell());
-                cell->get_dof_indices(local_dof_indices);
-
-                // delete a DoF index if it has not already been deleted
-                // (e.g., by visiting a neighboring cell, if it is on the
-                // boundary), and if we don't own it
-                for (types::global_dof_index &index : local_dof_indices)
-                  if ((index != numbers::invalid_dof_index) &&
-                      (!owned_dofs.is_element(index)))
-                    index = numbers::invalid_dof_index;
-
-                cell->set_dof_indices(local_dof_indices);
-              }
-        }
+        for (auto cell : dof_handler->active_cell_iterators())
+          if (cell->is_ghost())
+            {
+              DoFAccessorImplementation::Implementation::process_dof_indices(
+                *cell,
+                0,
+                cell->active_fe_index(),
+                DoFAccessorImplementation::Implementation::
+                  DoFIndexProcessor<dim, spacedim, false>(),
+                [&owned_dofs](auto &stored_index, auto) {
+                  // delete a DoF index if it has not already been
+                  // deleted (e.g., by visiting a neighboring cell, if
+                  // it is on the boundary), and if we don't own it
+                  if (stored_index != numbers::invalid_dof_index &&
+                      (!owned_dofs.is_element(stored_index)))
+                    stored_index = numbers::invalid_dof_index;
+                },
+                false);
+            }
 
 
         // renumber. Skip when there is nothing to do because we own no DoF.
@@ -4408,24 +4406,22 @@ namespace internal
 
         // delete all knowledge of DoF indices that are not locally
         // owned
-        {
-          std::vector<dealii::types::global_dof_index> local_dof_indices;
-
-          for (auto cell : dof_handler->cell_iterators_on_level(level))
-            if (cell->is_ghost_on_level())
-              {
-                local_dof_indices.resize(cell->get_fe().n_dofs_per_cell());
-                cell->get_mg_dof_indices(local_dof_indices);
-
-                for (types::global_dof_index &index : local_dof_indices)
-                  if ((index != numbers::invalid_dof_index) &&
-                      (!owned_dofs.is_element(index)))
-                    index = numbers::invalid_dof_index;
-
-                cell->set_mg_dof_indices(local_dof_indices);
-              }
-        }
-
+        for (auto cell : dof_handler->cell_iterators_on_level(level))
+          if (cell->is_ghost_on_level())
+            {
+              DoFAccessorImplementation::Implementation::process_dof_indices(
+                *cell,
+                0,
+                0,
+                DoFAccessorImplementation::Implementation::
+                  MGDoFIndexProcessor<dim, spacedim>(cell->level()),
+                [&owned_dofs](auto &stored_index, auto) {
+                  if ((stored_index != numbers::invalid_dof_index) &&
+                      (!owned_dofs.is_element(stored_index)))
+                    stored_index = numbers::invalid_dof_index;
+                },
+                true);
+            }
 
         // renumber. Skip when there is nothing to do because we own no DoF.
         if (level < triangulation->n_levels() && owned_dofs.n_elements() > 0)
index 9fdb50f4f9703dd258ba02f94120bf63c396e1b2..7681938603f1543fb62c1697da950d08ee8df370 100644 (file)
@@ -33,38 +33,6 @@ namespace internal
     {
       return (MemoryConsumption::memory_consumption(dofs));
     }
-
-
-
-    template <int dim>
-    template <int dh_dim, int spacedim>
-    void
-    DoFObjects<dim>::set_dof_index(
-      const dealii::DoFHandler<dh_dim, spacedim> &dof_handler,
-      const unsigned int                          obj_index,
-      const unsigned int                          fe_index,
-      const unsigned int                          local_index,
-      const types::global_dof_index               global_index)
-    {
-      (void)fe_index;
-      Assert(
-        (fe_index == dealii::DoFHandler<dh_dim, spacedim>::default_fe_index),
-        ExcMessage(
-          "Only the default FE index is allowed for non-hp-DoFHandler objects"));
-      Assert(
-        local_index < dof_handler.get_fe().template n_dofs_per_object<dim>(),
-        ExcIndexRange(local_index,
-                      0,
-                      dof_handler.get_fe().template n_dofs_per_object<dim>()));
-      Assert(obj_index *
-                   dof_handler.get_fe().template n_dofs_per_object<dim>() +
-                 local_index <
-               dofs.size(),
-             ExcInternalError());
-
-      dofs[obj_index * dof_handler.get_fe().template n_dofs_per_object<dim>() +
-           local_index] = global_index;
-    }
   } // namespace DoFHandlerImplementation
 } // namespace internal
 
index 9cf607ad1005311e7b2cd1b1dcd35b3cd80880b6..720f99f43e0ba2fa189b613d18339e52bbd66ff5 100644 (file)
@@ -19,55 +19,3 @@ for (deal_II_dimension : DIMENSIONS)
   {
     template class DoFObjects<deal_II_dimension>;
   }
-
-
-for (deal_II_dimension, structdim : DIMENSIONS)
-  {
-    template types::global_dof_index DoFObjects<structdim>::get_dof_index(
-      const dealii::DoFHandler<deal_II_dimension> &dof_handler,
-      const unsigned int                           obj_index,
-      const unsigned int                           fe_index,
-      const unsigned int                           local_index) const;
-
-    template void DoFObjects<structdim>::set_dof_index(
-      const dealii::DoFHandler<deal_II_dimension> &dof_handler,
-      const unsigned int                           obj_index,
-      const unsigned int                           fe_index,
-      const unsigned int                           local_index,
-      const types::global_dof_index                global_index);
-
-#if (deal_II_dimension < 3) && (structdim < 3)
-
-    template types::global_dof_index DoFObjects<structdim>::get_dof_index(
-      const dealii::DoFHandler<deal_II_dimension, deal_II_dimension + 1>
-        &                dof_handler,
-      const unsigned int obj_index,
-      const unsigned int fe_index,
-      const unsigned int local_index) const;
-
-    template void DoFObjects<structdim>::set_dof_index(
-      const dealii::DoFHandler<deal_II_dimension, deal_II_dimension + 1>
-        &                           dof_handler,
-      const unsigned int            obj_index,
-      const unsigned int            fe_index,
-      const unsigned int            local_index,
-      const types::global_dof_index global_index);
-#endif
-
-
-#if (deal_II_dimension == 3) && (structdim < 3)
-
-    template types::global_dof_index DoFObjects<structdim>::get_dof_index(
-      const dealii::DoFHandler<1, 3> &dof_handler,
-      const unsigned int              obj_index,
-      const unsigned int              fe_index,
-      const unsigned int              local_index) const;
-
-    template void DoFObjects<structdim>::set_dof_index(
-      const dealii::DoFHandler<1, 3> &dof_handler,
-      const unsigned int              obj_index,
-      const unsigned int              fe_index,
-      const unsigned int              local_index,
-      const types::global_dof_index   global_index);
-#endif
-  }

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.