]> https://gitweb.dealii.org/ - dealii.git/commitdiff
Check in DoFHandler::renumber_dofs() that DoFs are distributed 191/head
authorTimo Heister <timo.heister@gmail.com>
Thu, 9 Oct 2014 14:17:38 +0000 (10:17 -0400)
committerTimo Heister <timo.heister@gmail.com>
Fri, 10 Oct 2014 14:36:59 +0000 (10:36 -0400)
- This adds several Asserts to check that DoFs are distributed
before calling renumber_dofs() for active or level DoFs.
- Note that we now also require that you call distribute_dofs()
before distribute_mg_dofs().
- Updated documentation.

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

index 494257c39ab8c44771551a307c0772d95594db04..c38a675aa7c5be4f89889383b69a6870bc368090 100644 (file)
@@ -372,8 +372,11 @@ public:
   virtual void distribute_dofs (const FiniteElement<dim,spacedim> &fe);
 
   /**
-   * Distribute multigrid degrees of freedom similar
-   * to distribute_dofs() but on each level.
+   * Distribute level degrees of freedom on each level for geometric
+   * multigrid. The active DoFs need to be distributed using distribute_dofs()
+   * before calling this function and the @p fe needs to be identical to the
+   * finite element passed to distribute_dofs().
+   *
    * This replaces the functionality of the old MGDoFHandler.
    */
   virtual void distribute_mg_dofs (const FiniteElement<dim, spacedim> &fe);
index 94d240015dc7eafe85dbb6da20b53ac241eb7ddf..9f2c0380229d73d6f35ca19e380f61ad6ca12a5d 100644 (file)
@@ -1239,9 +1239,10 @@ void DoFHandler<dim,spacedim>::distribute_dofs (const FiniteElement<dim,spacedim
 template<int dim, int spacedim>
 void DoFHandler<dim, spacedim>::distribute_mg_dofs (const FiniteElement<dim, spacedim> &fe)
 {
+  Assert(levels.size()>0, ExcMessage("Distribute active DoFs using distribute_dofs() before calling distribute_mg_dofs()."));
+
   const FiniteElement<dim, spacedim> *old_fe = selected_fe;
-  Assert((old_fe==NULL) || (old_fe == &fe), ExcMessage("you are required to use the same FE for level and active DoFs!") );
-  selected_fe = &fe;
+  Assert(old_fe == &fe, ExcMessage("You are required to use the same FE for level and active DoFs!") );
 
   clear_mg_space();
 
@@ -1311,6 +1312,8 @@ template <int dim, int spacedim>
 void
 DoFHandler<dim,spacedim>::renumber_dofs (const std::vector<types::global_dof_index> &new_numbers)
 {
+  Assert(levels.size()>0, ExcMessage("You need to distribute DoFs before you can renumber them."));
+
   Assert (new_numbers.size() == n_locally_owned_dofs(),
           ExcRenumberingIncomplete());
 
@@ -1356,6 +1359,8 @@ template<>
 void DoFHandler<1>::renumber_dofs (const unsigned int level,
                                    const std::vector<types::global_dof_index> &new_numbers)
 {
+  Assert(mg_levels.size()>0 && levels.size()>0,
+         ExcMessage("You need to distribute active and level DoFs before you can renumber level DoFs."));
   Assert (new_numbers.size() == n_dofs(level), DoFHandler<1>::ExcRenumberingIncomplete());
 
   // note that we can not use cell iterators
@@ -1390,6 +1395,8 @@ template<>
 void DoFHandler<2>::renumber_dofs (const unsigned int  level,
                                    const std::vector<types::global_dof_index>  &new_numbers)
 {
+  Assert(mg_levels.size()>0 && levels.size()>0,
+         ExcMessage("You need to distribute active and level DoFs before you can renumber level DoFs."));
   Assert (new_numbers.size() == n_dofs(level),
           DoFHandler<2>::ExcRenumberingIncomplete());
 
@@ -1447,6 +1454,8 @@ template<>
 void DoFHandler<3>::renumber_dofs (const unsigned int  level,
                                    const std::vector<types::global_dof_index>  &new_numbers)
 {
+  Assert(mg_levels.size()>0 && levels.size()>0,
+         ExcMessage("You need to distribute active and level DoFs before you can renumber level DoFs."));
   Assert (new_numbers.size() == n_dofs(level),
           DoFHandler<3>::ExcRenumberingIncomplete());
 

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.