]> https://gitweb.dealii.org/ - dealii.git/commitdiff
Let DoFHandler policy classes return the NumberCache also in distribute_mg_dofs().
authorWolfgang Bangerth <bangerth@colostate.edu>
Fri, 23 Jun 2017 14:58:16 +0000 (08:58 -0600)
committerWolfgang Bangerth <bangerth@colostate.edu>
Fri, 23 Jun 2017 14:58:16 +0000 (08:58 -0600)
This follows the same pattern as a previous commit. However, it involves less
drama, with the exception of an innocent assertion that had to be removed
(see the last commit). Other than that, everything seems to work.

include/deal.II/dofs/dof_handler_policy.h
include/deal.II/dofs/number_cache.h
source/dofs/dof_handler.cc
source/dofs/dof_handler_policy.cc

index 4e725dfb6693e259a34b787f7def34ab9ae2f104..375ffd04b82f4aa839346bc7290532c900e1f487 100644 (file)
@@ -75,11 +75,12 @@ namespace internal
 
         /**
          * Distribute the multigrid dofs on each level of the DoFHandler
-         * associated with this policy object.
+         * associated with this policy object. Return a vector of number
+         * caches for all of the levels.
          */
         virtual
-        void
-        distribute_mg_dofs (std::vector<NumberCache> &number_caches) const = 0;
+        std::vector<NumberCache>
+        distribute_mg_dofs () const = 0;
 
         /**
          * Renumber degrees of freedom as specified by the first argument. The
@@ -118,8 +119,8 @@ namespace internal
 
         // documentation is inherited
         virtual
-        void
-        distribute_mg_dofs (std::vector<NumberCache> &number_caches) const;
+        std::vector<NumberCache>
+        distribute_mg_dofs () const;
 
         // documentation is inherited
         virtual
@@ -165,8 +166,8 @@ namespace internal
          * This function is not yet implemented.
          */
         virtual
-        void
-        distribute_mg_dofs (std::vector<NumberCache> &number_caches) const;
+        std::vector<NumberCache>
+        distribute_mg_dofs () const;
 
         /**
          * Renumber degrees of freedom as specified by the first argument.
@@ -206,8 +207,8 @@ namespace internal
 
         // documentation is inherited
         virtual
-        void
-        distribute_mg_dofs (std::vector<NumberCache> &number_caches) const;
+        std::vector<NumberCache>
+        distribute_mg_dofs () const;
 
         // documentation is inherited
         virtual
index 07615788e0e86e3f673ec2cc752033d84a2d05ad..16c44d7e2c86f49f2ff7d4083f729d6e8f1fc438 100644 (file)
@@ -39,12 +39,24 @@ namespace internal
        */
       NumberCache ();
 
+      /**
+       * Copy constructor. Simply copy all members of the referenced
+       * object to the current object.
+       */
+      NumberCache (const NumberCache &) = default;
+
       /**
        * Move constructor. Simply move all members of the referenced
        * object to the current object.
        */
       NumberCache (NumberCache &&) = default;
 
+      /**
+       * Copy operator. Simply copy all members of the referenced
+       * object to the current object.
+       */
+      NumberCache &operator= (const NumberCache &) = default;
+
       /**
        * Move assignment operator. Simply move all members of the referenced
        * object to the current object.
index 08f9011c18aa40f6d23059a89147cd4648e59acd..18be6f188273a3120d17e2dea06909d7e83abd2c 100644 (file)
@@ -1126,13 +1126,7 @@ void DoFHandler<dim, spacedim>::distribute_mg_dofs (const FiniteElement<dim, spa
   clear_mg_space();
 
   internal::DoFHandler::Implementation::reserve_space_mg (*this);
-  const parallel::distributed::Triangulation<dim,spacedim> *dist_tr = dynamic_cast<const parallel::distributed::Triangulation<dim,spacedim>*>(&*tria);
-  if (!dist_tr)
-    mg_number_cache.resize((*tria).n_levels());
-  else
-    mg_number_cache.resize(dist_tr->n_global_levels());
-
-  policy->distribute_mg_dofs (mg_number_cache);
+  mg_number_cache = policy->distribute_mg_dofs ();
 
   // initialize the block info object
   // only if this is a sequential
index 681056c601ef3d2cb1741350a72a49954eca11f3..86f87947354c3b799e53f571240fd8b4eba483cb 100644 (file)
@@ -957,29 +957,40 @@ namespace internal
 
 
       template <int dim, int spacedim>
-      void
+      std::vector<NumberCache>
       Sequential<dim,spacedim>::
-      distribute_mg_dofs (std::vector<NumberCache> &number_caches) const
+      distribute_mg_dofs () const
       {
         std::vector<bool> user_flags;
-
         dof_handler->get_triangulation().save_user_flags (user_flags);
+
         const_cast<dealii::Triangulation<dim, spacedim>&>(dof_handler->get_triangulation()).clear_user_flags ();
 
+        std::vector<NumberCache> number_caches;
+        number_caches.reserve (dof_handler->get_triangulation().n_levels());
         for (unsigned int level = 0; level < dof_handler->get_triangulation().n_levels(); ++level)
           {
-            types::global_dof_index next_free_dof = Implementation::distribute_dofs_on_level(0, numbers::invalid_subdomain_id,
-                                                    *dof_handler, level);
-
-            number_caches[level].n_global_dofs = next_free_dof;
-            number_caches[level].n_locally_owned_dofs = next_free_dof;
-            number_caches[level].locally_owned_dofs = complete_index_set(next_free_dof);
-            number_caches[level].locally_owned_dofs_per_processor.resize(1);
-            number_caches[level].locally_owned_dofs_per_processor[0] = complete_index_set(next_free_dof);
-            number_caches[level].n_locally_owned_dofs_per_processor.resize(1);
-            number_caches[level].n_locally_owned_dofs_per_processor[0] = next_free_dof;
+            const types::global_dof_index next_free_dof
+              = Implementation::distribute_dofs_on_level(0, numbers::invalid_subdomain_id,
+                                                         *dof_handler, level);
+
+            // set up the number cache of this level
+            NumberCache number_cache;
+            number_cache.n_global_dofs = next_free_dof;
+            number_cache.n_locally_owned_dofs = next_free_dof;
+            number_cache.locally_owned_dofs = complete_index_set(next_free_dof);
+            number_cache.locally_owned_dofs_per_processor.resize(1);
+            number_cache.locally_owned_dofs_per_processor[0] = complete_index_set(next_free_dof);
+            number_cache.n_locally_owned_dofs_per_processor.resize(1);
+            number_cache.n_locally_owned_dofs_per_processor[0] = next_free_dof;
+
+            // then push the current level's number cache onto the array
+            number_caches.emplace_back (number_cache);
           }
+
         const_cast<dealii::Triangulation<dim, spacedim>&>(dof_handler->get_triangulation()).load_user_flags (user_flags);
+
+        return number_caches;
       }
 
 
@@ -1270,14 +1281,17 @@ namespace internal
 
 
       template <int dim, int spacedim>
-      void
+      std::vector<NumberCache>
       ParallelShared<dim,spacedim>::
-      distribute_mg_dofs (std::vector<NumberCache> &number_caches) const
+      distribute_mg_dofs () const
       {
         // first, call the sequential function to distribute dofs
-        this->Sequential<dim,spacedim>::distribute_mg_dofs (number_caches);
+        std::vector<NumberCache> number_caches
+          = this->Sequential<dim,spacedim>::distribute_mg_dofs ();
         // now we need to update the number cache. This part is not yet implemented.
-        AssertThrow(false,ExcNotImplemented());
+        Assert(false,ExcNotImplemented());
+
+        return number_caches;
       }
 
 
@@ -2484,13 +2498,13 @@ namespace internal
 
 
       template <int dim, int spacedim>
-      void
+      std::vector<NumberCache>
       ParallelDistributed<dim, spacedim>::
-      distribute_mg_dofs (std::vector<NumberCache> &number_caches) const
+      distribute_mg_dofs () const
       {
 #ifndef DEAL_II_WITH_P4EST
-        (void)number_caches;
         Assert (false, ExcNotImplemented());
+        return std::vector<NumberCache>();
 #else
 
         parallel::distributed::Triangulation< dim, spacedim > *tr
@@ -2515,9 +2529,11 @@ namespace internal
         // Triangulation has fewer levels. we need to do this because
         // we need to communicate across all processors on all levels
         const unsigned int n_levels = tr->n_global_levels();
+        std::vector<NumberCache> number_caches;
+        number_caches.reserve(n_levels);
         for (unsigned int level = 0; level < n_levels; ++level)
           {
-            NumberCache &level_number_cache = number_caches[level];
+            NumberCache level_number_cache;
 
             //* 1. distribute on own subdomain
             const unsigned int n_initial_local_dofs =
@@ -2646,6 +2662,7 @@ namespace internal
                    == shift,
                    ExcInternalError());
 
+            number_caches.emplace_back (level_number_cache);
           }
 
 
@@ -2729,6 +2746,8 @@ namespace internal
         }
 #endif // DEBUG
 
+        return number_caches;
+
 #endif // DEAL_II_WITH_P4EST
       }
 

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.