/**
* The number of multilevel dofs on given level. Since hp::DoFHandler does
- * not support multilevel methods yet, this function returns
- * numbers::invalid_unsigned int independent of its argument.
+ * not support multilevel methods yet, this function throws an exception
+ * ExcNotImplemented() independent of its argument.
*/
types::global_dof_index n_dofs(const unsigned int level) const;
*/
const IndexSet &locally_owned_dofs() const;
-
/**
* Return a vector that stores the locally owned DoFs of each processor.
* If you are only interested in the number of elements each processor
const std::vector<types::global_dof_index> &
n_locally_owned_dofs_per_processor () const;
+ /**
+ * Return an IndexSet describing the set of locally owned DoFs used for
+ * the given multigrid level. Since hp::DoFHandler does not support
+ * multilevel methods yet, this function throws an exception
+ * ExcNotImplemented() independent of its argument.
+ */
+ const IndexSet &locally_owned_mg_dofs(const unsigned int level) const;
+
+ /**
+ * Return a vector that stores the locally owned level DoFs of each
+ * processor on the given level @p level. Since hp::DoFHandler does not
+ * support multilevel methods yet, this function throws an exception
+ * ExcNotImplemented() independent of its argument.
+ */
+ const std::vector<IndexSet> &
+ locally_owned_mg_dofs_per_processor (const unsigned int level) const;
+
/**
* Return a constant reference to the set of finite element objects that
* are used by this @p DoFHandler.
*/
dealii::internal::DoFHandler::NumberCache number_cache;
+ /**
+ * A structure that contains all sorts of numbers that characterize the
+ * degrees of freedom on multigrid levels. Since multigrid is not currently
+ * supported, this table is not filled with valid entries.
+ */
+ std::vector<dealii::internal::DoFHandler::NumberCache> mg_number_cache;
+
/**
* Array to store the indices for degrees of freedom located at vertices.
*
ar &vertex_dofs;
ar &vertex_dof_offsets;
ar &number_cache;
+ ar &mg_number_cache;
ar &levels;
ar &faces;
ar &has_children;
ar &vertex_dofs;
ar &vertex_dof_offsets;
ar &number_cache;
+ ar &mg_number_cache;
// boost::serialization can restore pointers just fine, but if the
// pointer object still points to something useful, that object is not
types::global_dof_index
DoFHandler<dim,spacedim>::n_dofs (const unsigned int) const
{
+ Assert(false, ExcNotImplemented());
return numbers::invalid_dof_index;
}
+ template <int dim, int spacedim>
+ const IndexSet &
+ DoFHandler<dim, spacedim>::locally_owned_mg_dofs (const unsigned int level) const
+ {
+ Assert(false, ExcNotImplemented());
+ (void)level;
+ Assert(level < this->get_triangulation().n_global_levels(),
+ ExcMessage("invalid level in locally_owned_mg_dofs"));
+ return mg_number_cache[0].locally_owned_dofs;
+ }
+
+
+ template <int dim, int spacedim>
+ const std::vector<IndexSet> &
+ DoFHandler<dim, spacedim>::locally_owned_mg_dofs_per_processor (const unsigned int level) const
+ {
+ Assert(false, ExcNotImplemented());
+ (void)level;
+ Assert(level < this->get_triangulation().n_global_levels(),
+ ExcMessage("invalid level in locally_owned_mg_dofs_per_processor"));
+ return mg_number_cache[0].locally_owned_dofs_per_processor;
+ }
+
+
+
template <int dim, int spacedim>
inline
const hp::FECollection<dim,spacedim> &
DoFTools::extract_locally_relevant_dofs<DoFHandler<deal_II_dimension-1,deal_II_dimension > >
(const DoFHandler<deal_II_dimension-1,deal_II_dimension > & dof_handler,
IndexSet & dof_set);
+ template
+ void
+ DoFTools::extract_locally_relevant_level_dofs<DoFHandler<deal_II_dimension-1,deal_II_dimension > >
+ (const DoFHandler<deal_II_dimension-1,deal_II_dimension > & dof_handler,
+ const unsigned int level,
+ IndexSet & dof_set);
+ template
+ void
+ DoFTools::extract_locally_relevant_level_dofs<hp::DoFHandler<deal_II_dimension-1,deal_II_dimension > >
+ (const hp::DoFHandler<deal_II_dimension-1,deal_II_dimension > & dof_handler,
+ const unsigned int level,
+ IndexSet & dof_set);
#endif
#if deal_II_dimension > 2
DoFTools::extract_locally_relevant_dofs<DoFHandler<deal_II_dimension-2,deal_II_dimension > >
(const DoFHandler<deal_II_dimension-2,deal_II_dimension > & dof_handler,
IndexSet & dof_set);
+ template
+ void
+ DoFTools::extract_locally_relevant_level_dofs<DoFHandler<deal_II_dimension-2,deal_II_dimension > >
+ (const DoFHandler<deal_II_dimension-2,deal_II_dimension > & dof_handler,
+ const unsigned int level,
+ IndexSet & dof_set);
+ template
+ void
+ DoFTools::extract_locally_relevant_level_dofs<hp::DoFHandler<deal_II_dimension-2,deal_II_dimension > >
+ (const hp::DoFHandler<deal_II_dimension-2,deal_II_dimension > & dof_handler,
+ const unsigned int level,
+ IndexSet & dof_set);
#endif
template
const unsigned int level,
IndexSet & dof_set);
+ template
+ void
+ DoFTools::extract_locally_relevant_level_dofs<hp::DoFHandler<deal_II_dimension> >
+ (const hp::DoFHandler<deal_II_dimension> & dof_handler,
+ const unsigned int level,
+ IndexSet & dof_set);
+
template
void
DoFTools::extract_constant_modes<DoFHandler<deal_II_dimension> >