]> https://gitweb.dealii.org/ - dealii-svn.git/commitdiff
Downstream numbering for DG
authorguido <guido@0785d39b-7218-0410-832d-ea1e28bc413d>
Thu, 13 Jul 2000 19:56:08 +0000 (19:56 +0000)
committerguido <guido@0785d39b-7218-0410-832d-ea1e28bc413d>
Thu, 13 Jul 2000 19:56:08 +0000 (19:56 +0000)
git-svn-id: https://svn.dealii.org/trunk@3158 0785d39b-7218-0410-832d-ea1e28bc413d

deal.II/deal.II/include/numerics/dof_renumbering.h
deal.II/deal.II/source/dofs/dof_renumbering.cc

index 6d3d61e495fa03780bb821aceaf9ae7887931b25..e4ad5aa13a2ec57d4b298a9b95d6469671377b6c 100644 (file)
  * degrees of freedom. This renumbering may only exchange whole blocks
  * and must not destroy the block structure.
  *
- * Given an ordered vector of cells, the function @p{cell_wise}
+ * Given an ordered vector of cells, the function @p{cell_wise_dg}
  * accomplishes this. Inside the cells, the previous ordering will be
  * preserved, so it may be useful to apply @component_wise} first.
  *
@@ -278,12 +278,44 @@ class DoFRenumbering
                                      * freedom have to be associated
                                      * with the interior of the cell.
                                      */
-    templateint dim>
+    template <int dim>
     static void
-    cell_wise (DoFHandler<dim>                     &dof_handler,
-              const vector<DoFCellAccessor<dim> > &cell_order);
+    cell_wise_dg (DoFHandler<dim>                     &dof_handler,
+                 const vector<typename DoFHandler<dim>::cell_iterator> &cell_order);
     
 
+                                    /**
+                                     * Cell-wise downstream numbering
+                                     * with respect to a constant
+                                     * flow direction.
+                                     *
+                                     * The cells are sorted such that
+                                     * the centers of higher numbers
+                                     * are further downstream with
+                                     * respect to the constant vector
+                                     * @p{direction} than the centers
+                                     * of lower numbers. Even if this
+                                     * yields a downstream numbering
+                                     * with respect to the flux on
+                                     * the edges for fairly general
+                                     * grids, this might not be
+                                     * guaranteed for all meshes.
+                                     *
+                                     * This function produces a
+                                     * downstream ordering of the
+                                     * mesh cells and calls
+                                     * @ref{cell_wise_dg}.
+                                     * Therefore, it only works with
+                                     * Discontinuous Galerkin Finite
+                                     * Elements, i.e. all degrees of
+                                     * freedom have to be associated
+                                     * with the interior of the cell.
+                                     */
+    template <int dim>
+    static void
+    downstream_dg (DoFHandler<dim>  &dof_handler,
+                  const Point<dim> &direction);
+
                                     /**
                                      * Sort those degrees of freedom
                                      * which are tagged with @p{true}
@@ -299,8 +331,8 @@ class DoFRenumbering
                                      */
     template <int dim>
     static void
-    sort_selected_dofs_back (const vector<bool> &selected_dofs,
-                            DoFHandler<dim>    &dof_handler);
+    sort_selected_dofs_back (DoFHandler<dim>    &dof_handler,
+                            const vector<bool> &selected_dofs);
 
                                     /**
                                      * Exception
index ae82208b61106729978a06a6bb2bd255042012cf..bb9ae94807a24f4bb51c494dc2a6c6d34c6ce29d 100644 (file)
@@ -501,8 +501,8 @@ void DoFRenumbering::component_wise (DoFHandler<dim>            &dof_handler,
 
 template <int dim>
 void
-DoFRenumbering::sort_selected_dofs_back (const vector<bool> &selected_dofs,
-                                        DoFHandler<dim>    &dof_handler)
+DoFRenumbering::sort_selected_dofs_back (DoFHandler<dim>    &dof_handler,
+                                        const vector<bool> &selected_dofs)
 {
   const unsigned int n_dofs = dof_handler.n_dofs();
   Assert (selected_dofs.size() == n_dofs,
@@ -538,8 +538,8 @@ DoFRenumbering::sort_selected_dofs_back (const vector<bool> &selected_dofs,
 
 template <int dim>
 void
-DoFRenumbering::cell_wise (DoFHandler<dim>& dof,
-                          const vector<DoFCellAccessor<dim> >& cells)
+DoFRenumbering::cell_wise_dg (DoFHandler<dim>& dof,
+                          const vector<typename DoFHandler<dim>::cell_iterator>& cells)
 {
   Assert(cells.size() == dof.get_tria().n_active_cells(),
         ExcDimensionMismatch(cells.size(),
@@ -565,18 +565,74 @@ DoFRenumbering::cell_wise (DoFHandler<dim>& dof,
 
   unsigned int global_index = 0;
   
-  typename vector<DoFCellAccessor<dim> >::const_iterator cell;
+  typename vector<typename DoFHandler<dim>::cell_iterator>::const_iterator cell;
 
   for(cell = cells.begin(); cell != cells.end(); ++cell)
     {
-      cell->get_dof_indices(cell_dofs);
+      (*cell)->get_dof_indices(cell_dofs);
       sort(cell_dofs.begin(), cell_dofs.end());
-      
+
       for (unsigned int i=0;i<n_cell_dofs;++i)
-       new_order[global_index++] = cell_dofs[i];
+       {
+         new_order[global_index++] = cell_dofs[i];
+       }
     }
   Assert(global_index == n_global_dofs, ExcRenumberingIncomplete());
-  dof.renumber_dofs(new_order);
+
+  vector<unsigned int> reverse(new_order.size());
+
+  for (unsigned int i=0;i<new_order.size(); ++i)
+    reverse[new_order[i]] = i;
+
+  dof.renumber_dofs(reverse);
+}
+
+/**
+ * Provide comparator for DoFCellAccessors
+ */
+
+template <int dim>
+struct CompCells
+{
+                                    /**
+                                     * Flow direction.
+                                     */
+    const Point<dim>& dir;
+                                    /**
+                                     * Constructor.
+                                     */
+    CompCells (const Point<dim>& dir) :
+                   dir(dir) 
+      {}
+                                    /**
+                                     * Return true if c1 < c2.
+                                     */    
+    bool operator () (const typename DoFHandler<dim>::cell_iterator& c1,
+                     const typename DoFHandler<dim>::cell_iterator&c2) const
+      {
+       Point<dim> diff = c2->center() - c1->center();
+       double s = diff*dir;
+       return (s>0);
+      }
+};
+
+    
+template <int dim>
+void
+DoFRenumbering::downstream_dg (DoFHandler<dim>& dof,
+                              const Point<dim>& direction)
+{
+  vector<typename DoFHandler<dim>::cell_iterator>
+    ordered_cells(dof.get_tria().n_active_cells());
+  CompCells<dim> comparator(direction);
+  
+  typename DoFHandler<dim>::active_cell_iterator begin = dof.begin_active();
+  typename DoFHandler<dim>::active_cell_iterator end = dof.end();
+  
+  copy (begin, end, ordered_cells.begin());
+  sort (ordered_cells.begin(), ordered_cells.end(), comparator);
+
+  cell_wise_dg(dof, ordered_cells);
 }
 
 
@@ -600,11 +656,13 @@ void DoFRenumbering::component_wise (DoFHandler<deal_II_dimension>&,
 
 
 template
-void DoFRenumbering::cell_wise (DoFHandler<deal_II_dimension>&,
-                               const vector<DoFCellAccessor<deal_II_dimension> >&);
-
+void DoFRenumbering::cell_wise_dg (DoFHandler<deal_II_dimension>&,
+                                  const vector<DoFHandler<deal_II_dimension>::cell_iterator>&);
 
+template
+void DoFRenumbering::downstream_dg (DoFHandler<deal_II_dimension>&,
+                                   const Point<deal_II_dimension>&);
 
 template
-void DoFRenumbering::sort_selected_dofs_back (const vector<bool> &,
-                                             DoFHandler<deal_II_dimension> &);
+void DoFRenumbering::sort_selected_dofs_back (DoFHandler<deal_II_dimension> &,
+                                             const vector<bool> &);

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.