]> https://gitweb.dealii.org/ - dealii-svn.git/commitdiff
fix DoFRenumbering::hierarchical for a distributed Triangulation
authorheister <heister@0785d39b-7218-0410-832d-ea1e28bc413d>
Thu, 25 Aug 2011 14:25:14 +0000 (14:25 +0000)
committerheister <heister@0785d39b-7218-0410-832d-ea1e28bc413d>
Thu, 25 Aug 2011 14:25:14 +0000 (14:25 +0000)
git-svn-id: https://svn.dealii.org/trunk@24194 0785d39b-7218-0410-832d-ea1e28bc413d

deal.II/include/deal.II/distributed/tria.h
deal.II/source/distributed/tria.cc
deal.II/source/dofs/dof_renumbering.cc

index 888afd1d6243f581bda94587b9c3a3693436af24..13c877498c3350e739e88a60d09adf93d88743b1 100644 (file)
@@ -495,6 +495,16 @@ namespace parallel
                                                                const CellStatus,
                                                                const void*)> & unpack_callback);
 
+                                        /**
+                                         * Returns a permutation vector for the order the coarse
+                                         * cells are handed of to p4est. For example the first
+                                         * element i in this vector denotes that the first cell
+                                         * in hierarchical ordering is the ith deal cell starting
+                                         * from begin(0).
+                                         */
+       const std::vector<unsigned int> &
+       get_p4est_tree_to_coarse_cell_permutation() const;
+
   private:
                                         /**
                                          * MPI communicator to be
@@ -699,6 +709,11 @@ namespace parallel
                                          */
        MPI_Comm get_communicator () const;
 
+       /**
+        * */
+               const std::vector<unsigned int> &
+       get_p4est_tree_to_coarse_cell_permutation() const;
+       
                                         /**
                                          * Return the subdomain id of
                                          * those cells that are owned
index 47191ee62563266c89d73bf0ce236a3b674adfc7..a548df84dea0dc0edfad7ade4a0c4b88fb2a13d4 100644 (file)
@@ -2803,6 +2803,15 @@ namespace parallel
 
 
 
+       template <int dim, int spacedim>
+       const std::vector<unsigned int> &
+       Triangulation<dim,spacedim>::get_p4est_tree_to_coarse_cell_permutation() const
+       {
+         return p4est_tree_to_coarse_cell_permutation;
+       }
+
+
+
     template <int dim, int spacedim>
     MPI_Comm
     Triangulation<dim,spacedim>::get_communicator () const
@@ -2943,6 +2952,13 @@ namespace parallel
       return MPI_COMM_WORLD;
     }
 
+       template <>
+       const std::vector<unsigned int> &
+       Triangulation<1,1>::get_p4est_tree_to_coarse_cell_permutation() const
+       {
+         static std::vector<unsigned int> a;
+         return a;
+       }
 
 
     template <>
index 7af3b620e2a30b2ec934372f809b01a57715eaab..9f306b1a53a5215352c7aa77c1cfeb59340adad3 100644 (file)
@@ -961,11 +961,42 @@ namespace DoFRenumbering
        unsigned int next_free = 0;
        const IndexSet locally_owned = dof_handler.locally_owned_dofs();
 
-    for(cell = dof_handler.begin(0); cell != dof_handler.end(0); ++cell)
-         next_free = compute_hierarchical_recursive<dim> (next_free,
-                                                                                                          renumbering,
-                                                                                                          cell,
-                                                                                                          locally_owned);
+    const parallel::distributed::Triangulation<dim> * tria
+       = dynamic_cast<const parallel::distributed::Triangulation<dim>*>
+          (&dof_handler.get_tria());
+
+       if (tria)
+      {
+#ifdef DEAL_II_USE_P4EST
+        //this is a distributed Triangulation. We need to traverse the coarse
+        //cells in the order p4est does
+        for (unsigned int c = 0; c < tria->n_cells (0); ++c)
+          {
+            unsigned int coarse_cell_index =
+              tria->get_p4est_tree_to_coarse_cell_permutation() [c];
+
+            const typename DoFHandler<dim>::cell_iterator
+            cell (tria, 0, coarse_cell_index, &dof_handler);
+
+            next_free = compute_hierarchical_recursive<dim> (next_free,
+                        renumbering,
+                        cell,
+                        locally_owned);
+          }
+#else
+        Assert (false, ExcNotImplemented());
+#endif
+      }
+    else
+      {
+        //this is not a distributed Triangulation. Traverse coarse cells in the
+        //normal order
+        for (cell = dof_handler.begin (0); cell != dof_handler.end (0); ++cell)
+          next_free = compute_hierarchical_recursive<dim> (next_free,
+                      renumbering,
+                      cell,
+                      locally_owned);
+      }
 
                                     // verify that the last numbered
                                     // degree of freedom is either
@@ -986,6 +1017,7 @@ namespace DoFRenumbering
                       numbers::invalid_unsigned_int)
            == renumbering.end(),
            ExcInternalError());
+       
        dof_handler.renumber_dofs(renumbering);
   }
 

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.