]> https://gitweb.dealii.org/ - dealii.git/commitdiff
Refactor DoFHandlerPolicy::Implementation::renumber_dofs().
authorWolfgang Bangerth <bangerth@colostate.edu>
Fri, 30 Jun 2017 15:32:39 +0000 (09:32 -0600)
committerWolfgang Bangerth <bangerth@colostate.edu>
Fri, 30 Jun 2017 18:12:44 +0000 (12:12 -0600)
The 1d, 2d, and 3d implementations of these functions had most of their code
duplicated. This can be done more elegantly by instead splitting the code into
vertices, cells, and faces. The main function then becomes dimension independent,
as are the functions dealing with vertices and cells, and only the face function
requires dimensional specialization.

source/dofs/dof_handler_policy.cc

index 2cd9303f1dea829aa397ae4cdad6dc7c3c40f979..e4fbc271b12fb2c9f95954a49143dce4b576b91e 100644 (file)
@@ -628,27 +628,20 @@ namespace internal
 
 
         /**
-         * Implementation of DoFHandler::renumber_dofs()
-         *
-         * If the second argument has any elements set, elements of
-         * the then the vector of new numbers do not relate to the old
-         * DoF number but instead to the index of the old DoF number
-         * within the set of locally owned DoFs.
+         * The part of the renumber_dofs() functionality that is dimension
+         * independent because it renumbers the DoF indices on vertices
+         * (which exist for all dimensions).
          *
-         * (The IndexSet argument is not used in 1d because we only need
-         * it for parallel meshes and 1d doesn't support that right now.)
+         * See renumber_dofs() for the meaning of the arguments.
          */
-        template <int spacedim>
+        template <int dim, int spacedim>
         static
         void
-        renumber_dofs (const std::vector<types::global_dof_index> &new_numbers,
-                       const IndexSet                             &indices,
-                       DoFHandler<1,spacedim>                     &dof_handler,
-                       const bool                                  check_validity)
+        renumber_vertex_dofs (const std::vector<types::global_dof_index> &new_numbers,
+                              const IndexSet                             &indices,
+                              DoFHandler<dim,spacedim>                   &dof_handler,
+                              const bool                                  check_validity)
         {
-          (void)indices;
-          Assert (indices == IndexSet(0), ExcNotImplemented());
-
           // we can not use cell iterators in this function since then
           // we would renumber the dofs on the interface of two cells
           // more than once. Anyway, this way it's not only more
@@ -659,7 +652,9 @@ namespace internal
                i=dof_handler.vertex_dofs.begin();
                i!=dof_handler.vertex_dofs.end(); ++i)
             if (*i != numbers::invalid_dof_index)
-              *i = new_numbers[*i];
+              *i = (indices.size() == 0)?
+                   (new_numbers[*i]) :
+                   (new_numbers[indices.index_within_set(*i)]);
             else if (check_validity)
               // if index is invalid_dof_index: check if this one
               // really is unused
@@ -668,13 +663,128 @@ namespace internal
                                    dof_handler.selected_fe->dofs_per_vertex)
                       == false,
                       ExcInternalError ());
+        }
+
 
+
+        /**
+         * The part of the renumber_dofs() functionality that is dimension
+         * independent because it renumbers the DoF indices on cell interiors
+         * (which exist for all dimensions).
+         *
+         * See renumber_dofs() for the meaning of the arguments.
+         */
+        template <int dim, int spacedim>
+        static
+        void
+        renumber_cell_dofs (const std::vector<types::global_dof_index> &new_numbers,
+                            const IndexSet                             &indices,
+                            DoFHandler<dim,spacedim>                   &dof_handler)
+        {
           for (unsigned int level=0; level<dof_handler.levels.size(); ++level)
             for (std::vector<types::global_dof_index>::iterator
                  i=dof_handler.levels[level]->dof_object.dofs.begin();
                  i!=dof_handler.levels[level]->dof_object.dofs.end(); ++i)
               if (*i != numbers::invalid_dof_index)
-                *i = new_numbers[*i];
+                *i = ((indices.size() == 0) ?
+                      new_numbers[*i] :
+                      new_numbers[indices.index_within_set(*i)]);
+        }
+
+
+
+        /**
+         * The part of the renumber_dofs() functionality that operates on faces.
+         * This part is dimension dependent and so needs to be implemented in
+         * three separate specializations of the function.
+         *
+         * See renumber_dofs() for the meaning of the arguments.
+         */
+        template <int spacedim>
+        static
+        void
+        renumber_face_dofs (const std::vector<types::global_dof_index> &/*new_numbers*/,
+                            const IndexSet                             &/*indices*/,
+                            DoFHandler<1,spacedim>                     &/*dof_handler*/)
+        {
+          // nothing to do in 1d since there are no separate faces
+        }
+
+
+
+        template <int spacedim>
+        static
+        void
+        renumber_face_dofs (const std::vector<types::global_dof_index> &new_numbers,
+                            const IndexSet                             &indices,
+                            DoFHandler<2,spacedim>                     &dof_handler)
+        {
+          // treat dofs on lines
+          for (std::vector<types::global_dof_index>::iterator
+               i=dof_handler.faces->lines.dofs.begin();
+               i!=dof_handler.faces->lines.dofs.end(); ++i)
+            if (*i != numbers::invalid_dof_index)
+              *i = ((indices.size() == 0) ?
+                    new_numbers[*i] :
+                    new_numbers[indices.index_within_set(*i)]);
+        }
+
+
+
+        template <int spacedim>
+        static
+        void
+        renumber_face_dofs (const std::vector<types::global_dof_index> &new_numbers,
+                            const IndexSet                             &indices,
+                            DoFHandler<3,spacedim>                     &dof_handler)
+        {
+          // treat dofs on lines
+          for (std::vector<types::global_dof_index>::iterator
+               i=dof_handler.faces->lines.dofs.begin();
+               i!=dof_handler.faces->lines.dofs.end(); ++i)
+            if (*i != numbers::invalid_dof_index)
+              *i = ((indices.size() == 0) ?
+                    new_numbers[*i] :
+                    new_numbers[indices.index_within_set(*i)]);
+
+          // treat dofs on quads
+          for (std::vector<types::global_dof_index>::iterator
+               i=dof_handler.faces->quads.dofs.begin();
+               i!=dof_handler.faces->quads.dofs.end(); ++i)
+            if (*i != numbers::invalid_dof_index)
+              *i = ((indices.size() == 0) ?
+                    new_numbers[*i] :
+                    new_numbers[indices.index_within_set(*i)]);
+        }
+
+
+
+
+        /**
+         * Implementation of DoFHandler::renumber_dofs()
+         *
+         * If the second argument has any elements set, elements of
+         * the then the vector of new numbers do not relate to the old
+         * DoF number but instead to the index of the old DoF number
+         * within the set of locally owned DoFs.
+         *
+         * (The IndexSet argument is not used in 1d because we only need
+         * it for parallel meshes and 1d doesn't support that right now.)
+         */
+        template <int dim, int spacedim>
+        static
+        void
+        renumber_dofs (const std::vector<types::global_dof_index> &new_numbers,
+                       const IndexSet                             &indices,
+                       DoFHandler<dim,spacedim>                   &dof_handler,
+                       const bool                                  check_validity)
+        {
+          if (dim == 1)
+            Assert (indices == IndexSet(0), ExcNotImplemented());
+
+          renumber_vertex_dofs (new_numbers, indices, dof_handler, check_validity);
+          renumber_face_dofs (new_numbers, indices, dof_handler);
+          renumber_cell_dofs (new_numbers, indices, dof_handler);
 
           // update the cache used for cell dof indices
           update_all_level_cell_dof_indices_caches(dof_handler);
@@ -682,6 +792,19 @@ namespace internal
 
 
 
+        template <int dim, int spacedim>
+        static
+        void
+        renumber_dofs (const std::vector<types::global_dof_index> &/*new_numbers*/,
+                       const IndexSet                             &/*indices*/,
+                       hp::DoFHandler<dim,spacedim>               &/*dof_handler*/,
+                       const bool                                  /*check_validity*/)
+        {
+          Assert (false, ExcNotImplemented());
+        }
+
+
+
         template <int spacedim>
         static
         void
@@ -732,61 +855,6 @@ namespace internal
 
 
 
-        template <int spacedim>
-        static
-        void
-        renumber_dofs (const std::vector<types::global_dof_index> &new_numbers,
-                       const IndexSet                             &indices,
-                       DoFHandler<2,spacedim>                     &dof_handler,
-                       const bool                                  check_validity)
-        {
-          // we can not use cell iterators in this function since then
-          // we would renumber the dofs on the interface of two cells
-          // more than once. Anyway, this way it's not only more
-          // correct but also faster; note, however, that dof numbers
-          // may be invalid_dof_index, namely when the appropriate
-          // vertex/line/etc is unused
-          for (std::vector<types::global_dof_index>::iterator
-               i=dof_handler.vertex_dofs.begin();
-               i!=dof_handler.vertex_dofs.end(); ++i)
-            if (*i != numbers::invalid_dof_index)
-              *i = (indices.size() == 0)?
-                   (new_numbers[*i]) :
-                   (new_numbers[indices.index_within_set(*i)]);
-            else if (check_validity)
-              // if index is invalid_dof_index: check if this one
-              // really is unused
-              Assert (dof_handler.get_triangulation()
-                      .vertex_used((i-dof_handler.vertex_dofs.begin()) /
-                                   dof_handler.selected_fe->dofs_per_vertex)
-                      == false,
-                      ExcInternalError ());
-
-          for (std::vector<types::global_dof_index>::iterator
-               i=dof_handler.faces->lines.dofs.begin();
-               i!=dof_handler.faces->lines.dofs.end(); ++i)
-            if (*i != numbers::invalid_dof_index)
-              *i = ((indices.size() == 0) ?
-                    new_numbers[*i] :
-                    new_numbers[indices.index_within_set(*i)]);
-
-          for (unsigned int level=0; level<dof_handler.levels.size(); ++level)
-            {
-              for (std::vector<types::global_dof_index>::iterator
-                   i=dof_handler.levels[level]->dof_object.dofs.begin();
-                   i!=dof_handler.levels[level]->dof_object.dofs.end(); ++i)
-                if (*i != numbers::invalid_dof_index)
-                  *i = ((indices.size() == 0) ?
-                        new_numbers[*i] :
-                        new_numbers[indices.index_within_set(*i)]);
-            }
-
-          // update the cache used for cell dof indices
-          update_all_level_cell_dof_indices_caches (dof_handler);
-        }
-
-
-
         template <int spacedim>
         static
         void
@@ -871,68 +939,6 @@ namespace internal
 
 
 
-        template <int spacedim>
-        static
-        void
-        renumber_dofs (const std::vector<types::global_dof_index> &new_numbers,
-                       const IndexSet                             &indices,
-                       DoFHandler<3,spacedim>                     &dof_handler,
-                       const bool                                  check_validity)
-        {
-          // we can not use cell iterators in this function since then
-          // we would renumber the dofs on the interface of two cells
-          // more than once. Anyway, this way it's not only more
-          // correct but also faster; note, however, that dof numbers
-          // may be invalid_dof_index, namely when the appropriate
-          // vertex/line/etc is unused
-          for (std::vector<types::global_dof_index>::iterator
-               i=dof_handler.vertex_dofs.begin();
-               i!=dof_handler.vertex_dofs.end(); ++i)
-            if (*i != numbers::invalid_dof_index)
-              *i = ((indices.size() == 0) ?
-                    new_numbers[*i] :
-                    new_numbers[indices.index_within_set(*i)]);
-            else if (check_validity)
-              // if index is invalid_dof_index: check if this one
-              // really is unused
-              Assert (dof_handler.get_triangulation()
-                      .vertex_used((i-dof_handler.vertex_dofs.begin()) /
-                                   dof_handler.selected_fe->dofs_per_vertex)
-                      == false,
-                      ExcInternalError ());
-
-          for (std::vector<types::global_dof_index>::iterator
-               i=dof_handler.faces->lines.dofs.begin();
-               i!=dof_handler.faces->lines.dofs.end(); ++i)
-            if (*i != numbers::invalid_dof_index)
-              *i = ((indices.size() == 0) ?
-                    new_numbers[*i] :
-                    new_numbers[indices.index_within_set(*i)]);
-          for (std::vector<types::global_dof_index>::iterator
-               i=dof_handler.faces->quads.dofs.begin();
-               i!=dof_handler.faces->quads.dofs.end(); ++i)
-            if (*i != numbers::invalid_dof_index)
-              *i = ((indices.size() == 0) ?
-                    new_numbers[*i] :
-                    new_numbers[indices.index_within_set(*i)]);
-
-          for (unsigned int level=0; level<dof_handler.levels.size(); ++level)
-            {
-              for (std::vector<types::global_dof_index>::iterator
-                   i=dof_handler.levels[level]->dof_object.dofs.begin();
-                   i!=dof_handler.levels[level]->dof_object.dofs.end(); ++i)
-                if (*i != numbers::invalid_dof_index)
-                  *i = ((indices.size() == 0) ?
-                        new_numbers[*i] :
-                        new_numbers[indices.index_within_set(*i)]);
-            }
-
-          // update the cache used for cell dof indices
-          update_all_level_cell_dof_indices_caches (dof_handler);
-        }
-
-
-
         template <int spacedim>
         static
         void
@@ -1044,19 +1050,6 @@ namespace internal
 
 
 
-        template <int dim, int spacedim>
-        static
-        void
-        renumber_dofs (const std::vector<types::global_dof_index> &/*new_numbers*/,
-                       const IndexSet                             &/*indices*/,
-                       hp::DoFHandler<dim,spacedim>               &/*dof_handler*/,
-                       const bool                                  /*check_validity*/)
-        {
-          Assert (false, ExcNotImplemented());
-        }
-
-
-
         template <int dim, int spacedim>
         static
         void

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.