Reference documentation for deal.II version 9.0.0
|
#include <deal.II/fe/component_mask.h>
Public Member Functions | |
ComponentMask ()=default | |
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 ::ExceptionBase & | ExcNoComponentSelected () |
Private Attributes | |
std::vector< bool > | component_mask |
Friends | |
std::ostream & | operator<< (std::ostream &out, const ComponentMask &mask) |
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.
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):
The result is a component mask that, in 2d, would have values [false, false, true]
. Similarly, using
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.
|
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.
|
inline |
Initialize an object of this type with a set of selected components specified by the argument.
component_mask | A 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. |
Definition at line 255 of file component_mask.h.
|
inline |
Initialize the component mask with a number of elements that are either all true or false.
n_components | The number of elements of this mask |
initializer | The value each of these elements is supposed to have: either true or false. |
Definition at line 262 of file component_mask.h.
|
inline |
Set a particular entry in the mask to a value.
Definition at line 279 of file component_mask.h.
|
inline |
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.
Definition at line 271 of file component_mask.h.
|
inline |
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.
component_index | The 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. |
Definition at line 288 of file component_mask.h.
|
inline |
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.
Definition at line 306 of file component_mask.h.
|
inline |
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.
Definition at line 316 of file component_mask.h.
|
inline |
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.
Definition at line 342 of file component_mask.h.
|
inline |
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.
Definition at line 364 of file component_mask.h.
|
inline |
Return a component mask that contains the union of the components selected by the current object and the one passed as an argument.
Definition at line 373 of file component_mask.h.
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.
|
inline |
Return whether this object and the argument are identical.
Definition at line 421 of file component_mask.h.
|
inline |
Return whether this object and the argument are not identical.
Definition at line 429 of file component_mask.h.
std::size_t ComponentMask::memory_consumption | ( | ) | const |
Determine an estimate for the memory consumption (in bytes) of this object.
Definition at line 47 of file component_mask.cc.
|
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]
.
out | The stream to write to. |
mask | The mask to write. |
Definition at line 24 of file component_mask.cc.
|
private |
The actual component mask.
Definition at line 228 of file component_mask.h.