Reference documentation for deal.II version GIT relicensing-216-gcda843ca70 2024-03-28 12:20:02+00:00
\(\newcommand{\dealvcentcolon}{\mathrel{\mathop{:}}}\) \(\newcommand{\dealcoloneq}{\dealvcentcolon\mathrel{\mkern-1.2mu}=}\) \(\newcommand{\jump}[1]{\left[\!\left[ #1 \right]\!\right]}\) \(\newcommand{\average}[1]{\left\{\!\left\{ #1 \right\}\!\right\}}\)
Loading...
Searching...
No Matches
Public Member Functions | Static Public Member Functions | Private Attributes | Friends | List of all members

#include <deal.II/fe/component_mask.h>

Public Member Functions

 ComponentMask ()=default
 
template<typename = void>
 ComponentMask (const std::vector< bool > &block_mask)
 
 ComponentMask (const std::vector< bool > &component_mask)
 
 ComponentMask (const unsigned int n_components, const bool initializer)
 
void set (const unsigned int index, const bool value)
 
unsigned int size () const
 
bool operator[] (const unsigned int component_index) const
 
bool represents_n_components (const unsigned int n) const
 
unsigned int n_selected_components (const unsigned int overall_number_of_components=numbers::invalid_unsigned_int) const
 
unsigned int first_selected_component (const unsigned int overall_number_of_components=numbers::invalid_unsigned_int) const
 
bool represents_the_all_selected_mask () const
 
ComponentMask operator| (const ComponentMask &mask) const
 
ComponentMask operator& (const ComponentMask &mask) const
 
bool operator== (const ComponentMask &mask) const
 
bool operator!= (const ComponentMask &mask) const
 
std::size_t memory_consumption () const
 

Static Public Member Functions

static ::ExceptionBaseExcNoComponentSelected ()
 

Private Attributes

std::vector< boolcomponent_mask
 

Friends

std::ostream & operator<< (std::ostream &out, const ComponentMask &mask)
 

Detailed Description

This class represents a mask that can be used to select individual vector components of a finite element (see also this glossary entry). It will typically have as many elements as the finite element has vector components, and one can use operator[] to query whether a particular component has been selected.

Note
A "mask" represents a data structure with true and false entries that is generally used to enable or disable an operation for a particular vector component. By this definition, disabled vector components still exist – they are simply not touched. As a consequence, when you apply a component mask for interpolating boundary values (to choose just one example) of a problem with \(C\) vector components, the input argument that describes the boundary values will still have to provide \(C\) components even if the mask says that we only want to interpolate a subset of these components onto the finite element space. In other words, a component mask does not represent a reduction operation; it represents a selection.

Objects of this kind are used in many places where one wants to restrict operations to a certain subset of components, e.g. in DoFTools::make_zero_boundary_values() or VectorTools::interpolate_boundary_values(). These objects can either be created by hand, or, simpler, by asking the finite element to generate a component mask from certain selected components using code such as this where we create a mask that only denotes the velocity components of a Stokes element (see Handling vector valued problems):

// Q2 element for the velocities, Q1 element for the pressure
FESystem<dim> stokes_fe (FE_Q<dim>(2), dim,
FE_Q<dim>(1), 1);
ComponentMask pressure_mask = stokes_fe.component_mask (pressure);
std::vector< bool > component_mask
Definition fe_q.h:550

The result is a component mask that, in 2d, would have values [false, false, true]. Similarly, using

ComponentMask velocity_mask = stokes_fe.component_mask (velocities);

would result in a mask [true, true, false] in 2d. Of course, in 3d, the result would be [true, true, true, false].

Definition at line 80 of file component_mask.h.

Constructor & Destructor Documentation

◆ ComponentMask() [1/4]

ComponentMask::ComponentMask ( )
default

Initialize a component mask. The default is that a component mask represents a set of components that are all selected, i.e., calling this constructor results in a component mask that always returns true whenever asked whether a component is selected.

◆ ComponentMask() [2/4]

template<typename = void>
ComponentMask::ComponentMask ( const std::vector< bool > &  block_mask)

Deprecated constructor allowing implicit conversion.

◆ ComponentMask() [3/4]

ComponentMask::ComponentMask ( const std::vector< bool > &  component_mask)
explicit

Initialize an object of this type with a set of selected components specified by the argument.

Parameters
component_maskA vector of true/false entries that determine which components of a finite element are selected. If the length of the given vector is zero, then this interpreted as the case where every component is selected.

◆ ComponentMask() [4/4]

ComponentMask::ComponentMask ( const unsigned int  n_components,
const bool  initializer 
)

Initialize the component mask with a number of elements that are either all true or false.

Parameters
n_componentsThe number of elements of this mask
initializerThe value each of these elements is supposed to have: either true or false.

Member Function Documentation

◆ set()

void ComponentMask::set ( const unsigned int  index,
const bool  value 
)

Set a particular entry in the mask to a value.

◆ size()

unsigned int ComponentMask::size ( ) const

If this component mask has been initialized with a mask of size greater than zero, then return the size of the mask represented by this object. On the other hand, if this mask has been initialized as an empty object that represents a mask that is true for every element (i.e., if this object would return true when calling represents_the_all_selected_mask()) then return zero since no definite size is known.

◆ operator[]()

bool ComponentMask::operator[] ( const unsigned int  component_index) const

Return whether a particular component is selected by this mask. If this mask represents the case of an object that selects all components (e.g. if it is created using the default constructor or is converted from an empty vector of type bool) then this function returns true regardless of the given argument.

Parameters
component_indexThe index for which the function should return whether the component is selected. If this object represents a mask in which all components are always selected then any index is allowed here. Otherwise, the given index needs to be between zero and the number of components that this mask represents.

◆ represents_n_components()

bool ComponentMask::represents_n_components ( const unsigned int  n) const

Return whether this component mask represents a mask with exactly n components. This is true if either it was initialized with a vector with exactly n entries of type bool (in this case, n must equal the result of size()) or if it was initialized with an empty vector (or using the default constructor) in which case it can represent a mask with an arbitrary number of components and will always say that a component is selected.

◆ n_selected_components()

unsigned int ComponentMask::n_selected_components ( const unsigned int  overall_number_of_components = numbers::invalid_unsigned_int) const

Return the number of components that are selected by this mask.

Since empty component masks represent a component mask that would return true for every component, this function may not know the true size of the component mask and it therefore requires an argument that denotes the overall number of components.

If the object has been initialized with a non-empty mask (i.e., if the size() function returns something greater than zero, or equivalently if represents_the_all_selected_mask() returns false) then the argument can be omitted and the result of size() is taken.

◆ first_selected_component()

unsigned int ComponentMask::first_selected_component ( const unsigned int  overall_number_of_components = numbers::invalid_unsigned_int) const

Return the index of the first selected component. The argument is there for the same reason it exists with the n_selected_components() function.

The function throws an exception if no component is selected at all.

◆ represents_the_all_selected_mask()

bool ComponentMask::represents_the_all_selected_mask ( ) const

Return true if this mask represents a default constructed mask that corresponds to one in which all components are selected. If true, then the size() function will return zero.

◆ operator|()

ComponentMask ComponentMask::operator| ( const ComponentMask mask) const

Return a component mask that contains the union of the components selected by the current object and the one passed as an argument.

◆ operator&()

ComponentMask ComponentMask::operator& ( const ComponentMask mask) const

Return a component mask that has only those elements set that are set both in the current object as well as the one passed as an argument.

◆ operator==()

bool ComponentMask::operator== ( const ComponentMask mask) const

Return whether this object and the argument are identical.

◆ operator!=()

bool ComponentMask::operator!= ( const ComponentMask mask) const

Return whether this object and the argument are not identical.

◆ memory_consumption()

std::size_t ComponentMask::memory_consumption ( ) const

Determine an estimate for the memory consumption (in bytes) of this object.

Definition at line 46 of file component_mask.cc.

Friends And Related Symbol Documentation

◆ operator<<

std::ostream & operator<< ( std::ostream &  out,
const ComponentMask mask 
)
friend

Write a component mask to an output stream. If the component mask represents one where all components are selected without specifying a particular size of the mask, then it writes the string [all components selected] to the stream. Otherwise, it prints the component mask in a form like [true,true,true,false].

Parameters
outThe stream to write to.
maskThe mask to write.
Returns
A reference to the first argument.

Definition at line 23 of file component_mask.cc.

Member Data Documentation

◆ component_mask

std::vector<bool> ComponentMask::component_mask
private

The actual component mask.

Definition at line 244 of file component_mask.h.


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