* may be difficult, however, and in many cases will not justify the effort.
*
*
- * @sect2{Componentwise numbering}
+ * @sect2{Component-wise numbering}
*
* For finite elements composed of several base elements using the @p{FESystem}
* class, or for elements which provide several components themselves, it
* algorithm and afterwards renumbering component-wise. This will bring out the
* matrix structure and additionally have a good numbering within each block.
*
+ * @sect2{Cell-wise numbering for Discontinuous Galerkin FEM}
+ *
+ * One advantage of DGFEM is the fact, that it yields invertible
+ * blocks on the diagonal of the global matrix. these blocks are in
+ * fact the cell matrices and matrix elements outside these blocks are
+ * due to fluxes.
+ *
+ * Still, it may be necessary to apply a downstream numbering of the
+ * 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}
+ * accomplishes this. Inside the cells, the previous ordering will be
+ * preserved, so it may be useful to apply @component_wise} first.
*
* @sect2{Multigrid DoF numbering}
*
* to the actual function declarations to get more information on this.
*
*
- * @author Wolfgang Bangerth, 1998, 1999
+ * @author Wolfgang Bangerth, Guido Kanschat, 1998, 1999, 2000
*/
class DoFRenumbering
{
component_wise (DoFHandler<dim> &dof_handler,
const vector<unsigned int> &component_order = vector<unsigned int>());
+ /**
+ * Cell-wise renumbering for DG
+ * elements. This function takes
+ * the ordered set of cells in
+ * @p{cell_order}, and makes sure
+ * that all degrees of freedom in
+ * a cell with higher index are
+ * behind all degrees of freedom
+ * of a cell with lower
+ * index. The order inside a cell
+ * bloock will be the same as
+ * before this renumbering.
+ *
+ * This function 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
+ cell_wise (DoFHandler<dim> &dof_handler,
+ const vector<DoFCellAccessor<dim> > &cell_order);
+
+
/**
* Sort those degrees of freedom
* which are tagged with @p{true}
* Exception
*/
DeclException0 (ExcInvalidComponentOrder);
+
/**
- * Exception
+ * The function is only
+ * implemented for Discontinuous
+ * Galerkin Finite elements.
*/
- DeclException2 (ExcInvalidArraySize,
- int, int,
- << "The array has size " << arg1 << " but should have "
- << arg2 << ".");
+ DeclException0 (ExcNotDGFEM);
};
#include <lac/sparsity_pattern.h>
#include <dofs/dof_accessor.h>
#include <grid/tria_iterator.h>
+#include <grid/tria.h>
#include <dofs/dof_handler.h>
#include <dofs/dof_constraints.h>
#include <dofs/dof_tools.h>
{
const unsigned int n_dofs = dof_handler.n_dofs();
Assert (selected_dofs.size() == n_dofs,
- ExcInvalidArraySize (selected_dofs.size(), n_dofs));
+ ExcDimensionMismatch (selected_dofs.size(), n_dofs));
// re-sort the dofs according to
// their selection state
};
+template <int dim>
+void
+DoFRenumbering::cell_wise (DoFHandler<dim>& dof,
+ const vector<DoFCellAccessor<dim> >& cells)
+{
+ Assert(cells.size() == dof.get_tria().n_active_cells(),
+ ExcDimensionMismatch(cells.size(),
+ dof.get_tria().n_active_cells()));
+ switch (dim)
+ {
+ case 3:
+ Assert(dof.get_fe().n_dofs_per_quad()==0,
+ ExcNotDGFEM());
+ case 2:
+ Assert(dof.get_fe().n_dofs_per_line()==0,
+ ExcNotDGFEM());
+ default:
+ Assert(dof.get_fe().n_dofs_per_vertex()==0,
+ ExcNotDGFEM());
+ }
+
+ unsigned int n_global_dofs = dof.n_dofs();
+ unsigned int n_cell_dofs = dof.get_fe().n_dofs_per_cell();
+
+ vector<unsigned int> new_order(n_global_dofs);
+ vector<unsigned int> cell_dofs(n_cell_dofs);
+
+ unsigned int global_index = 0;
+
+ typename vector<DoFCellAccessor<dim> >::const_iterator cell;
+
+ for(cell = cells.begin(); cell != cells.end(); ++cell)
+ {
+ 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];
+ }
+ Assert(global_index == n_global_dofs, ExcRenumberingIncomplete());
+ dof.renumber_dofs(new_order);
+}
// explicit instantiations
template
-void DoFRenumbering::Cuthill_McKee (DoFHandler<deal_II_dimension> &dof_handler,
- const bool reversed_numbering,
- const bool use_constraints,
- const vector<unsigned int> &starting_indices);
+void DoFRenumbering::Cuthill_McKee (DoFHandler<deal_II_dimension>&,
+ const bool,
+ const bool,
+ const vector<unsigned int>&);
+
+template
+void DoFRenumbering::Cuthill_McKee (MGDoFHandler<deal_II_dimension>&,
+ const unsigned int,
+ const bool,
+ const vector<unsigned int>&);
template
-void DoFRenumbering::Cuthill_McKee (MGDoFHandler<deal_II_dimension> &dof_handler,
- const unsigned int level,
- const bool reversed_numbering,
- const vector<unsigned int> &starting_indices);
+void DoFRenumbering::component_wise (DoFHandler<deal_II_dimension>&,
+ const vector<unsigned int>&);
+
template
-void DoFRenumbering::component_wise (DoFHandler<deal_II_dimension> &dof_handler,
- const vector<unsigned int> &component_order_arg);
+void DoFRenumbering::cell_wise (DoFHandler<deal_II_dimension>&,
+ const vector<DoFCellAccessor<deal_II_dimension> >&);
+
template
-void DoFRenumbering::sort_selected_dofs_back (const vector<bool> &selected_dofs,
- DoFHandler<deal_II_dimension> &dof_handler);
+void DoFRenumbering::sort_selected_dofs_back (const vector<bool> &,
+ DoFHandler<deal_II_dimension> &);