Reference documentation for deal.II version 9.1.1
\(\newcommand{\dealcoloneq}{\mathrel{\vcenter{:}}=}\)
Public Member Functions | Private Attributes | List of all members
MGConstrainedDoFs Class Reference

#include <deal.II/multigrid/mg_constrained_dofs.h>

Inheritance diagram for MGConstrainedDoFs:
[legend]

Public Member Functions

template<int dim, int spacedim>
void initialize (const DoFHandler< dim, spacedim > &dof)
 
template<int dim, int spacedim>
void initialize (const DoFHandler< dim, spacedim > &dof, const std::map< types::boundary_id, const Function< spacedim > *> &function_map, const ComponentMask &component_mask=ComponentMask())
 
template<int dim, int spacedim>
void make_zero_boundary_constraints (const DoFHandler< dim, spacedim > &dof, const std::set< types::boundary_id > &boundary_ids, const ComponentMask &component_mask=ComponentMask())
 
template<int dim, int spacedim>
void make_no_normal_flux_constraints (const DoFHandler< dim, spacedim > &dof, const types::boundary_id bid, const unsigned int first_vector_component)
 
void clear ()
 
bool is_boundary_index (const unsigned int level, const types::global_dof_index index) const
 
bool at_refinement_edge (const unsigned int level, const types::global_dof_index index) const
 
bool is_interface_matrix_entry (const unsigned int level, const types::global_dof_index i, const types::global_dof_index j) const
 
const IndexSetget_boundary_indices (const unsigned int level) const
 
const IndexSetget_refinement_edge_indices (unsigned int level) const
 
bool have_boundary_indices () const
 
const AffineConstraints< double > & get_level_constraints (const unsigned int level) const
 
const AffineConstraints< double > & get_level_constraint_matrix (const unsigned int level) const
 
- Public Member Functions inherited from Subscriptor
 Subscriptor ()
 
 Subscriptor (const Subscriptor &)
 
 Subscriptor (Subscriptor &&) noexcept
 
virtual ~Subscriptor ()
 
Subscriptoroperator= (const Subscriptor &)
 
Subscriptoroperator= (Subscriptor &&) noexcept
 
void subscribe (std::atomic< bool > *const validity, const std::string &identifier="") const
 
void unsubscribe (std::atomic< bool > *const validity, const std::string &identifier="") const
 
unsigned int n_subscriptions () const
 
template<typename StreamType >
void list_subscribers (StreamType &stream) const
 
void list_subscribers () const
 
template<class Archive >
void serialize (Archive &ar, const unsigned int version)
 

Private Attributes

std::vector< IndexSetboundary_indices
 
std::vector< IndexSetrefinement_edge_indices
 
std::vector< AffineConstraints< double > > level_constraints
 

Additional Inherited Members

- Static Public Member Functions inherited from Subscriptor
static ::ExceptionBaseExcInUse (int arg1, std::string arg2, std::string arg3)
 
static ::ExceptionBaseExcNoSubscriber (std::string arg1, std::string arg2)
 

Detailed Description

Collection of boundary constraints and refinement edge constraints for level vectors.

Definition at line 42 of file mg_constrained_dofs.h.

Member Function Documentation

◆ initialize() [1/2]

template<int dim, int spacedim>
void MGConstrainedDoFs::initialize ( const DoFHandler< dim, spacedim > &  dof)
inline

Fill the internal data structures with hanging node constraints extracted from the dof handler object. Works with natural boundary conditions only. There exists a sister function setting up boundary constraints as well.

This function ensures that on every level, degrees of freedom at interior edges of a refinement level are treated corrected but leaves degrees of freedom at the boundary of the domain untouched assuming that no Dirichlet boundary conditions for them exist.

Furthermore, this call sets up an AffineConstraints object on each level that contains possible periodicity constraints in case those have been added to the underlying triangulation. The AffineConstraints object can be queried by get_level_constraints(level). Note that the current implementation of periodicity constraints in this class does not support rotation matrices in the periodicity definition, i.e., the respective argument in the GridTools::collect_periodic_faces() may not be different from the identity matrix.

Definition at line 218 of file mg_constrained_dofs.h.

◆ initialize() [2/2]

template<int dim, int spacedim>
void MGConstrainedDoFs::initialize ( const DoFHandler< dim, spacedim > &  dof,
const std::map< types::boundary_id, const Function< spacedim > *> &  function_map,
const ComponentMask component_mask = ComponentMask() 
)
inline

Fill the internal data structures with values extracted from the dof handler object and apply the boundary values provided.

This function internally calls the initialize() function above and the constrains degrees on the external boundary of the domain by calling MGTools::make_boundary_list() with the given second and third argument.

Deprecated:
Use initialize() followed by make_zero_boundary_constraints() instead

Definition at line 303 of file mg_constrained_dofs.h.

◆ make_zero_boundary_constraints()

template<int dim, int spacedim>
void MGConstrainedDoFs::make_zero_boundary_constraints ( const DoFHandler< dim, spacedim > &  dof,
const std::set< types::boundary_id > &  boundary_ids,
const ComponentMask component_mask = ComponentMask() 
)
inline

Fill the internal data structures with information about Dirichlet boundary dofs.

The initialize() function has to be called before to set hanging node constraints.

This function can be called multiple times to allow considering different sets of boundary_ids for different components.

Definition at line 325 of file mg_constrained_dofs.h.

◆ make_no_normal_flux_constraints()

template<int dim, int spacedim>
void MGConstrainedDoFs::make_no_normal_flux_constraints ( const DoFHandler< dim, spacedim > &  dof,
const types::boundary_id  bid,
const unsigned int  first_vector_component 
)
inline

Fill the internal data structures with information about no normal flux boundary dofs.

This function is limited to meshes whose no normal flux boundaries have faces which are normal to the x-, y-, or z-axis. Also, for a specific boundary id, all faces must be facing in the same direction, i.e., a boundary normal to the x-axis must have a different boundary id than a boundary normal to the y- or z-axis and so on. If the mesh was produced, for example, using the GridGenerator::hyper_cube() function, setting colorize=true during mesh generation and calling make_no_normal_flux_constraints() for each no normal flux boundary is sufficient.

Definition at line 346 of file mg_constrained_dofs.h.

◆ clear()

void MGConstrainedDoFs::clear ( )
inline

Reset the data structures.

Definition at line 400 of file mg_constrained_dofs.h.

◆ is_boundary_index()

bool MGConstrainedDoFs::is_boundary_index ( const unsigned int  level,
const types::global_dof_index  index 
) const
inline

Determine whether a dof index is subject to a boundary constraint.

Definition at line 408 of file mg_constrained_dofs.h.

◆ at_refinement_edge()

bool MGConstrainedDoFs::at_refinement_edge ( const unsigned int  level,
const types::global_dof_index  index 
) const
inline

Determine whether a dof index is at the refinement edge.

Definition at line 419 of file mg_constrained_dofs.h.

◆ is_interface_matrix_entry()

bool MGConstrainedDoFs::is_interface_matrix_entry ( const unsigned int  level,
const types::global_dof_index  i,
const types::global_dof_index  j 
) const
inline

Determine whether the (i,j) entry of the interface matrix on a given level should be set. This is taken in terms of dof i, that is, return true if i is at a refinement edge, j is not, and both are not on the external boundary.

Definition at line 428 of file mg_constrained_dofs.h.

◆ get_boundary_indices()

const IndexSet & MGConstrainedDoFs::get_boundary_indices ( const unsigned int  level) const
inline

Return the indices of level dofs on the given level that are subject to Dirichlet boundary conditions (as set by the function_map parameter in initialize()). The indices are restricted to the set of locally relevant level dofs.

Definition at line 445 of file mg_constrained_dofs.h.

◆ get_refinement_edge_indices()

const IndexSet & MGConstrainedDoFs::get_refinement_edge_indices ( unsigned int  level) const
inline

Return the indices of dofs on the given level that lie on an refinement edge (dofs on faces to neighbors that are coarser).

Definition at line 454 of file mg_constrained_dofs.h.

◆ have_boundary_indices()

bool MGConstrainedDoFs::have_boundary_indices ( ) const
inline

Return if Dirichlet boundary indices are set in initialize().

Definition at line 463 of file mg_constrained_dofs.h.

◆ get_level_constraints()

const AffineConstraints< double > & MGConstrainedDoFs::get_level_constraints ( const unsigned int  level) const
inline

Return the AffineConstraints object for a given level, containing periodicity constraints (if enabled on the triangulation).

Definition at line 471 of file mg_constrained_dofs.h.

◆ get_level_constraint_matrix()

const AffineConstraints< double > & MGConstrainedDoFs::get_level_constraint_matrix ( const unsigned int  level) const
inline

Return the AffineConstraints object for a given level, containing periodicity constraints (if enabled on the triangulation).

Deprecated:
Use get_level_constraints instead, which has a more descriptive name.

Definition at line 480 of file mg_constrained_dofs.h.

Member Data Documentation

◆ boundary_indices

std::vector<IndexSet> MGConstrainedDoFs::boundary_indices
private

The indices of boundary dofs for each level.

Definition at line 200 of file mg_constrained_dofs.h.

◆ refinement_edge_indices

std::vector<IndexSet> MGConstrainedDoFs::refinement_edge_indices
private

The degrees of freedom on a given level that live on the refinement edge between the level and cells on a coarser level.

Definition at line 206 of file mg_constrained_dofs.h.

◆ level_constraints

std::vector<AffineConstraints<double> > MGConstrainedDoFs::level_constraints
private

Constraint matrices containing information regarding potential periodic boundary conditions for each level .

Definition at line 212 of file mg_constrained_dofs.h.


The documentation for this class was generated from the following file: