16 #ifndef dealii_mg_constrained_dofs_h 17 #define dealii_mg_constrained_dofs_h 19 #include <deal.II/base/config.h> 20 #include <deal.II/base/subscriptor.h> 22 #include <deal.II/multigrid/mg_tools.h> 27 DEAL_II_NAMESPACE_OPEN
29 template <
int dim,
int spacedim>
class DoFHandler;
30 template <
int dim,
typename Number>
struct FunctionMap;
43 typedef std::vector<std::set<types::global_dof_index> >::size_type size_dof;
63 template <
int dim,
int spacedim>
76 template <
int dim,
int spacedim>
92 template <
int dim,
int spacedim>
94 const std::set<types::boundary_id> &boundary_ids,
177 template <
int dim,
int spacedim>
191 for (
unsigned int l=0; l<nlevels; ++l)
201 for (; cell!=endc; ++cell)
204 for (
unsigned int f = 0;
205 f < GeometryInfo<dim>::faces_per_cell;
207 if (cell->has_periodic_neighbor(f))
209 if (cell->is_locally_owned_on_level())
212 ExcMessage (
"Periodic neighbor of a locally owned cell must either be owned or ghost."));
222 const unsigned int dofs_per_face = cell->face(f)->
get_fe(0).dofs_per_face;
223 std::vector<types::global_dof_index> dofs_1(dofs_per_face);
224 std::vector<types::global_dof_index> dofs_2(dofs_per_face);
226 cell->periodic_neighbor(f)->face(cell->periodic_neighbor_face_no(f))->get_mg_dof_indices(l, dofs_1, 0);
227 cell->face(f)->get_mg_dof_indices(l, dofs_2, 0);
232 for (
unsigned int i=0; i<dofs_per_face; ++i)
253 template <
int dim,
int spacedim>
275 template <
int dim,
int spacedim>
279 const std::set<types::boundary_id> &boundary_ids,
333 const IndexSet &interface_dofs_on_level
389 DEAL_II_NAMESPACE_CLOSE
cell_iterator begin(const unsigned int level=0) const
bool is_boundary_index(const unsigned int level, const types::global_dof_index index) const
cell_iterator end() const
#define AssertIndexRange(index, range)
std::vector< IndexSet > boundary_indices
std::vector< ConstraintMatrix > level_constraints
bool is_interface_matrix_entry(const unsigned int level, const types::global_dof_index i, const types::global_dof_index j) const
static ::ExceptionBase & ExcMessage(std::string arg1)
const IndexSet & get_boundary_indices(const unsigned int level) const
ActiveSelector::cell_iterator cell_iterator
unsigned int global_dof_index
#define Assert(cond, exc)
const FiniteElement< dim, spacedim > & get_fe(const unsigned int index=0) const
const ConstraintMatrix & get_level_constraint_matrix(const unsigned int level) const
types::global_dof_index n_dofs() const
std::map< types::boundary_id, const Function< dim, Number > * > type
const IndexSet & get_refinement_edge_indices(unsigned int level) const
void make_zero_boundary_constraints(const DoFHandler< dim, spacedim > &dof, const std::set< types::boundary_id > &boundary_ids, const ComponentMask &component_mask=ComponentMask())
std::vector< IndexSet > refinement_edge_indices
const types::subdomain_id artificial_subdomain_id
bool at_refinement_edge(const unsigned int level, const types::global_dof_index index) const
const Triangulation< dim, spacedim > & get_triangulation() const
bool is_element(const size_type index) const
bool have_boundary_indices() const
static ::ExceptionBase & ExcInternalError()
void initialize(const DoFHandler< dim, spacedim > &dof)