#include <lac/forward_declarations.h>
#include <grid/forward_declarations.h>
+//TODO: Consider incorporating these functions in DoFTools
/**
* This is a collection of functions operating on, and manipulating
* for more information.
*
* All member functions are static, so there is no need to create an
- * object of class #DoFTools#.
+ * object of class #MGDoFTools#.
*
* @author Wolfgang Bangerth and others, 1999
*/
/**
* Write the sparsity structure
* of the matrix belonging to the
- * specified #level# including
- * constrained degrees of freedom
- * into the matrix structure. The
- * sparsity pattern does not
- * include entries introduced by
- * the elimination of constrained
- * nodes. The sparsity pattern
- * is not compressed, since if
- * you want to call
- * #ConstraintMatrix::condense(1)#
- * afterwards, new entries have
- * to be added. However, if you
- * don't want to call
- * #ConstraintMatrix::condense(1)#,
+ * specified #level#. The sparsity pattern
+ * is not compressed, so before
+ * creating the actual matrix
* you have to compress the
* matrix yourself, using
* #SparseMatrixStruct::compress()#.
+ *
+ * There is no need to consider
+ * hanging nodes here, since only
+ * one level is considered.
*/
template <int dim>
static void
make_sparsity_pattern (const MGDoFHandler<dim> &dof_handler,
- const unsigned int level,
- SparsityPattern &sparsity);
+ SparsityPattern &sparsity,
+ const unsigned int level);
+
+ /**
+ * Make a sparsity pattern including fluxes
+ * of discontinuous Galerkin methods.
+ * @see make_sparsity_pattern
+ * $see DoFTools
+ */
+ template <int dim>
+ static void
+ make_flux_sparsity_pattern (const MGDoFHandler<dim> &dof_handler,
+ SparsityPattern &sparsity,
+ const unsigned int level);
+
+
};
template <int dim>
-void MGDoFTools::make_sparsity_pattern (const MGDoFHandler<dim> &mg_dof_handler,
- const unsigned int level,
- SparsityPattern &sparsity)
+void MGDoFTools::make_sparsity_pattern (const MGDoFHandler<dim> &dof,
+ SparsityPattern &sparsity,
+ const unsigned int level)
{
- Assert (sparsity.n_rows() == mg_dof_handler.n_dofs(level),
- ExcDimensionMismatch (sparsity.n_rows(), mg_dof_handler.n_dofs(level)));
- Assert (sparsity.n_cols() == mg_dof_handler.n_dofs(level),
- ExcDimensionMismatch (sparsity.n_cols(), mg_dof_handler.n_dofs(level)));
+ const unsigned int n_dofs = dof.n_dofs(level);
- const unsigned int dofs_per_cell = mg_dof_handler.get_fe().dofs_per_cell;
+ Assert (sparsity.n_rows() == n_dofs,
+ ExcDimensionMismatch (sparsity.n_rows(), n_dofs));
+ Assert (sparsity.n_cols() == n_dofs,
+ ExcDimensionMismatch (sparsity.n_cols(), n_dofs));
+
+ const unsigned int dofs_per_cell = dof.get_fe().dofs_per_cell;
vector<unsigned int> dofs_on_this_cell(dofs_per_cell);
- MGDoFHandler<dim>::cell_iterator cell = mg_dof_handler.begin(level),
- endc = mg_dof_handler.end(level);
+ MGDoFHandler<dim>::cell_iterator cell = dof.begin(level),
+ endc = dof.end(level);
for (; cell!=endc; ++cell)
{
cell->get_mg_dof_indices (dofs_on_this_cell);
for (unsigned int j=0; j<dofs_per_cell; ++j)
sparsity.add (dofs_on_this_cell[i],
dofs_on_this_cell[j]);
- };
-};
+ }
+}
+
+template<int dim>
+void
+MGDoFTools::make_flux_sparsity_pattern (const MGDoFHandler<dim> &dof,
+ SparsityPattern &sparsity,
+ const unsigned int level)
+{
+ const unsigned int n_dofs = dof.n_dofs(level);
+
+ Assert (sparsity.n_rows() == n_dofs,
+ ExcDimensionMismatch (sparsity.n_rows(), n_dofs));
+ Assert (sparsity.n_cols() == n_dofs,
+ ExcDimensionMismatch (sparsity.n_cols(), n_dofs));
+
+ const unsigned int dofs_per_cell = dof.get_fe().dofs_per_cell;
+ vector<unsigned int> dofs_on_this_cell(dofs_per_cell);
+ vector<unsigned int> dofs_on_other_cell(dofs_per_cell);
+ DoFHandler<dim>::cell_iterator cell = dof.begin(level),
+ endc = dof.end(level);
+ for (; cell!=endc; ++cell)
+ {
+ cell->get_dof_indices (dofs_on_this_cell);
+ // make sparsity pattern for this cell
+ for (unsigned int i=0; i<dofs_per_cell; ++i)
+ for (unsigned int j=0; j<dofs_per_cell; ++j)
+ sparsity.add (dofs_on_this_cell[i],
+ dofs_on_this_cell[j]);
+ // Loop over all interior neighbors
+ for (unsigned int face = 0;
+ face < GeometryInfo<dim>::faces_per_cell;
+ ++face)
+ {
+ if ( (! cell->at_boundary(face)) && (cell->neighbor_level(face) == level) )
+ {
+ DoFHandler<dim>::cell_iterator neighbor = cell->neighbor(face);
+ neighbor->get_dof_indices (dofs_on_other_cell);
+ // only add one direction
+ // The other is taken care of
+ // by neighbor.
+ for (unsigned int i=0; i<dofs_per_cell; ++i)
+ {
+ for (unsigned int j=0; j<dofs_per_cell; ++j)
+ {
+ sparsity.add (dofs_on_this_cell[i],
+ dofs_on_other_cell[j]);
+ }
+ }
+ }
+ }
+ }
+}
// explicit instantiations
template void MGDoFTools::make_sparsity_pattern (const MGDoFHandler<deal_II_dimension> &,
- const unsigned int,
- SparsityPattern &);
+ SparsityPattern &,
+ const unsigned int);
+template void MGDoFTools::make_flux_sparsity_pattern (const MGDoFHandler<deal_II_dimension> &,
+ SparsityPattern &,
+ const unsigned int);