]> https://gitweb.dealii.org/ - dealii.git/commitdiff
Simplify code by using new NumberCache constructors.
authorWolfgang Bangerth <bangerth@colostate.edu>
Sat, 24 Jun 2017 08:22:40 +0000 (02:22 -0600)
committerWolfgang Bangerth <bangerth@colostate.edu>
Mon, 26 Jun 2017 05:02:47 +0000 (23:02 -0600)
source/dofs/dof_handler_policy.cc

index 52cbae372ea0805b59a3bcda5ad02d8af1329456..1e9551337132a8b3fc74a0cb42b7b2bc71e6a6f0 100644 (file)
@@ -930,27 +930,8 @@ namespace internal
                                            numbers::invalid_subdomain_id,
                                            *dof_handler);
 
-        // now set the elements of the
-        // number cache appropriately
-        NumberCache number_cache;
-        number_cache.n_global_dofs        = n_dofs;
-        number_cache.n_locally_owned_dofs = number_cache.n_global_dofs;
-
-        number_cache.locally_owned_dofs
-          = IndexSet (number_cache.n_global_dofs);
-        number_cache.locally_owned_dofs.add_range (0,
-                                                   number_cache.n_global_dofs);
-        number_cache.locally_owned_dofs.compress();
-
-        number_cache.n_locally_owned_dofs_per_processor
-          = std::vector<types::global_dof_index> (1,
-                                                  number_cache.n_global_dofs);
-
-        number_cache.locally_owned_dofs_per_processor
-          = std::vector<IndexSet> (1,
-                                   number_cache.locally_owned_dofs);
-
-        return number_cache;
+        // return a sequential, complete index set
+        return NumberCache (n_dofs);
       }
 
 
@@ -969,22 +950,13 @@ namespace internal
         number_caches.reserve (dof_handler->get_triangulation().n_levels());
         for (unsigned int level = 0; level < dof_handler->get_triangulation().n_levels(); ++level)
           {
-            const types::global_dof_index next_free_dof
+            // first distribute dofs on this level
+            const types::global_dof_index n_level_dofs
               = Implementation::distribute_dofs_on_level(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);
+            // then add a complete, sequential index set
+            number_caches.emplace_back (NumberCache(n_level_dofs));
           }
 
         const_cast<dealii::Triangulation<dim, spacedim>&>(dof_handler->get_triangulation()).load_user_flags (user_flags);
@@ -1002,26 +974,8 @@ namespace internal
         Implementation::renumber_dofs (new_numbers, IndexSet(0),
                                        *dof_handler, true);
 
-        // in the sequential case, the number cache should not have changed but we
-        // have to set the elements of the structure appropriately anyway
-        NumberCache number_cache;
-        number_cache.n_global_dofs        = dof_handler->n_dofs();
-        number_cache.n_locally_owned_dofs = number_cache.n_global_dofs;
-
-        number_cache.locally_owned_dofs
-          = IndexSet (number_cache.n_global_dofs);
-        number_cache.locally_owned_dofs.add_range (0,
-                                                   number_cache.n_global_dofs);
-        number_cache.locally_owned_dofs.compress();
-
-        number_cache.n_locally_owned_dofs_per_processor
-          = std::vector<types::global_dof_index> (1,
-                                                  number_cache.n_global_dofs);
-
-        number_cache.locally_owned_dofs_per_processor
-          = std::vector<IndexSet> (1,
-                                   number_cache.locally_owned_dofs);
-        return number_cache;
+        // return a sequential, complete index set
+        return NumberCache (new_numbers.size());
       }
 
 
@@ -1146,7 +1100,7 @@ namespace internal
 
         // first let the sequential algorithm do its magic. it is going to
         // enumerate DoFs on all cells, regardless of owner
-        NumberCache number_cache
+        NumberCache sequential_number_cache
           = this->Sequential<dim,spacedim>::distribute_dofs ();
 
         // then re-enumerate them based on their subdomain association.
@@ -1158,7 +1112,7 @@ namespace internal
         // function has to deal with other kinds of triangulations as
         // well, whereas we here know what kind of triangulation
         // we have and can simplify the code accordingly
-        const unsigned int n_dofs = number_cache.n_global_dofs;
+        const unsigned int n_dofs = sequential_number_cache.n_global_dofs;
         std::vector<types::global_dof_index> new_dof_indices (n_dofs,
                                                               numbers::invalid_dof_index);
         {
@@ -1194,7 +1148,9 @@ namespace internal
         // version of the function because we do things on all
         // cells and all cells have their subdomain ids and DoFs
         // correctly set
-        number_cache = this->Sequential<dim, spacedim>::renumber_dofs (new_dof_indices);
+        //
+        // the return value is not of interest to us here
+        this->Sequential<dim, spacedim>::renumber_dofs (new_dof_indices);
 
 
         // update the number cache. for this, we first have to find the subdomain
@@ -1213,11 +1169,9 @@ namespace internal
           Assert (subdomain_association[i] >= subdomain_association[i-1],
                   ExcInternalError());
 
-        {
-          number_cache.locally_owned_dofs_per_processor.clear();
-          number_cache.locally_owned_dofs_per_processor.resize (n_procs,
+        std::vector<IndexSet> locally_owned_dofs_per_processor (n_procs,
                                                                 IndexSet(n_dofs));
-
+        {
           // we know that the set of subdomain indices is contiguous from
           // the assertion above; find the start and end index for each
           // processor, taking into account that sometimes a processor
@@ -1239,7 +1193,7 @@ namespace internal
               if (end_index > start_index)
                 {
                   const unsigned int subdomain_owner = subdomain_association[start_index];
-                  number_cache.locally_owned_dofs_per_processor[subdomain_owner]
+                  locally_owned_dofs_per_processor[subdomain_owner]
                   .add_range (start_index, end_index);
                 }
 
@@ -1248,20 +1202,6 @@ namespace internal
             }
         }
 
-        // from this, calculate all of the other number cache fields
-        number_cache.locally_owned_dofs
-          = number_cache.locally_owned_dofs_per_processor[this->dof_handler->get_triangulation().locally_owned_subdomain()];
-        number_cache.n_locally_owned_dofs_per_processor
-        .resize (number_cache.locally_owned_dofs_per_processor.size());
-        for (unsigned int i=0; i<number_cache.n_locally_owned_dofs_per_processor.size(); ++i)
-          number_cache.n_locally_owned_dofs_per_processor[i]
-            = number_cache.locally_owned_dofs_per_processor[i].n_elements();
-        number_cache.n_locally_owned_dofs
-          = number_cache.n_locally_owned_dofs_per_processor[this->dof_handler->get_triangulation().locally_owned_subdomain()];
-
-        number_cache.n_global_dofs = n_dofs;
-
-
         // finally, restore current subdomain ids
         if (tr->with_artificial_cells())
           {
@@ -1273,7 +1213,10 @@ namespace internal
               cell->set_subdomain_id(saved_subdomain_ids[index]);
           }
 
-        return number_cache;
+        // return a NumberCache object made up from the sets of locally
+        // owned DoFs
+        return NumberCache(locally_owned_dofs_per_processor,
+                           this->dof_handler->get_triangulation().locally_owned_subdomain ());
       }
 
 
@@ -1284,7 +1227,7 @@ namespace internal
       distribute_mg_dofs () const
       {
         // first, call the sequential function to distribute dofs
-        std::vector<NumberCache> number_caches
+        const 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.
         Assert(false,ExcNotImplemented());
@@ -1411,22 +1354,13 @@ namespace internal
                    ExcInternalError());
           }
 
-        NumberCache number_cache
-          = this->Sequential<dim, spacedim>::renumber_dofs (global_gathered_numbers);
-        // correct number_cache:
-        number_cache.locally_owned_dofs_per_processor =
-          DoFTools::locally_owned_dofs_per_subdomain (*this->dof_handler);
-        number_cache.locally_owned_dofs =
-          number_cache.locally_owned_dofs_per_processor[this->dof_handler->get_triangulation().locally_owned_subdomain ()];
-        // sequential renumbering returns a vector of size 1 here,
-        // correct this:
-        number_cache.n_locally_owned_dofs_per_processor.resize(number_cache.locally_owned_dofs_per_processor.size());
-        for (unsigned int i = 0;
-             i < number_cache.n_locally_owned_dofs_per_processor.size (); i++)
-          number_cache.n_locally_owned_dofs_per_processor[i] = number_cache.locally_owned_dofs_per_processor[i].n_elements ();
-
-        number_cache.n_locally_owned_dofs =
-          number_cache.n_locally_owned_dofs_per_processor[this->dof_handler->get_triangulation().locally_owned_subdomain ()];
+        // let the sequential algorithm do its magic; ignore the
+        // return type, but reconstruct the number cache based on
+        // which DoFs each process owns
+        this->Sequential<dim, spacedim>::renumber_dofs (global_gathered_numbers);
+
+        const NumberCache number_cache (DoFTools::locally_owned_dofs_per_subdomain (*this->dof_handler),
+                                        this->dof_handler->get_triangulation().locally_owned_subdomain ());
 
         // restore artificial cells
         cell = tr->begin_active();
@@ -2575,7 +2509,7 @@ namespace internal
                     }
               }
 
-
+// TODO: make this code simpler with the new constructors of NumberCache
             // make indices consecutive
             level_number_cache.n_locally_owned_dofs = 0;
             for (std::vector<dealii::types::global_dof_index>::iterator it=renumbering.begin();

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.