]> https://gitweb.dealii.org/ - dealii.git/commitdiff
Limit setting and getting active_fe_indices to cells where we can do that.
authorWolfgang Bangerth <bangerth@colostate.edu>
Sat, 22 Jul 2017 00:54:29 +0000 (18:54 -0600)
committerWolfgang Bangerth <bangerth@colostate.edu>
Sun, 13 Aug 2017 00:17:10 +0000 (18:17 -0600)
include/deal.II/dofs/dof_accessor.templates.h
source/hp/dof_handler.cc

index 371ee3d994bc21f0ea52778e121475fad885bc12..b7c55d4eeee9abc27283b8bccc3ab45323ae2711 100644 (file)
@@ -3073,8 +3073,7 @@ namespace internal
       unsigned int
       active_fe_index (const DoFCellAccessor<DoFHandler<dim,spacedim>, level_dof_access> &)
       {
-        // ::DoFHandler only supports a
-        // single active fe with index
+        // ::DoFHandler only supports a single active fe with index
         // zero
         return 0;
       }
@@ -3106,8 +3105,7 @@ namespace internal
                            const unsigned int                                i)
       {
         (void)i;
-        // ::DoFHandler only supports a
-        // single active fe with index
+        // ::DoFHandler only supports a single active fe with index
         // zero
         typedef dealii::DoFAccessor<dim,DoFHandler<dim,spacedim>, level_dof_access> BaseClass;
         Assert (i == 0, typename BaseClass::ExcInvalidObject());
@@ -3757,6 +3755,14 @@ DoFCellAccessor<DoFHandlerType,level_dof_access>::active_fe_index () const
                       "children because no degrees of freedom are assigned "
                       "to this cell and, consequently, no finite element "
                       "is associated with it."));
+  Assert ((dynamic_cast<const dealii::DoFHandler<DoFHandlerType::dimension,DoFHandlerType::space_dimension>*>
+           (this->dof_handler) != nullptr)
+          ||
+          (this->is_locally_owned() || this->is_ghost()),
+          ExcMessage ("You can only query active_fe_index information on cells "
+                      "that are either locally owned or (after distributing "
+                      "degrees of freedom) are ghost cells."));
+
   return dealii::internal::DoFCellAccessor::Implementation::active_fe_index (*this);
 }
 
@@ -3774,6 +3780,17 @@ DoFCellAccessor<DoFHandlerType,level_dof_access>::set_active_fe_index (const uns
           ExcMessage ("You can not set the active_fe_index on a cell that has "
                       "children because no degrees of freedom will be assigned "
                       "to this cell."));
+
+  if (i != DoFHandlerType::default_fe_index)
+    Assert ((dynamic_cast<const dealii::DoFHandler<DoFHandlerType::dimension,DoFHandlerType::space_dimension>*>
+             (this->dof_handler) != nullptr)
+            ||
+            this->is_locally_owned(),
+            ExcMessage ("You can only set active_fe_index information on cells "
+                        "that are locally owned. On ghost cells, this information "
+                        "will automatically be propagated from the owning process "
+                        "of that cell, and there is no information at all on "
+                        "artificial cells."));
   dealii::internal::DoFCellAccessor::Implementation::set_active_fe_index (*this, i);
 }
 
index a37f1798fb62867d8a5409d9fab1a0ec9469661e..fef78f22f1f8eb60606a984f6bb54fae78567bdd 100644 (file)
@@ -1237,7 +1237,8 @@ namespace hp
     active_cell_iterator cell=begin_active(),
                          endc=end();
     for (unsigned int i=0; cell!=endc; ++cell, ++i)
-      cell->set_active_fe_index(active_fe_indices[i]);
+      if (cell->is_locally_owned())
+        cell->set_active_fe_index(active_fe_indices[i]);
   }
 
 
@@ -1253,7 +1254,7 @@ namespace hp
     active_cell_iterator cell=begin_active(),
                          endc=end();
     for (unsigned int i=0; cell!=endc; ++cell, ++i)
-      active_fe_indices[i]=cell->active_fe_index();
+      active_fe_indices[i] = cell->active_fe_index();
   }
 
 
@@ -1603,9 +1604,14 @@ namespace hp
                 // cell->active_fe_index() since that function is not
                 // allowed for inactive cells, but we can access this
                 // information from the DoFLevels directly
+                //
+                // we don't have to set the active_fe_index for ghost
+                // cells -- these will be exchanged automatically
+                // upon distribute_dofs()
                 for (unsigned int i = 0; i < cell->n_children(); ++i)
-                  cell->child (i)->set_active_fe_index
-                  (levels[cell->level()]->active_fe_index (cell->index()));
+                  if (cell->child (i)->is_locally_owned())
+                    cell->child (i)->set_active_fe_index
+                    (levels[cell->level()]->active_fe_index (cell->index()));
               }
           }
       }

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.