Reference documentation for deal.II version GIT relicensing-214-g6e74dec06b 2024-03-27 18:10:01+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 | Private Attributes | Friends | List of all members

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

Public Member Functions

 BlockMask ()=default
 
template<typename = void>
 BlockMask (const std::vector< bool > &block_mask)
 
 BlockMask (const std::vector< bool > &block_mask)
 
 BlockMask (const unsigned int n_blocks, const bool initializer)
 
unsigned int size () const
 
bool operator[] (const unsigned int block_index) const
 
bool represents_n_blocks (const unsigned int n) const
 
unsigned int n_selected_blocks (const unsigned int overall_number_of_blocks=numbers::invalid_unsigned_int) const
 
unsigned int first_selected_block (const unsigned int overall_number_of_blocks=numbers::invalid_unsigned_int) const
 
bool represents_the_all_selected_mask () const
 
BlockMask operator| (const BlockMask &mask) const
 
BlockMask operator& (const BlockMask &mask) const
 
bool operator== (const BlockMask &mask) const
 
bool operator!= (const BlockMask &mask) const
 
std::size_t memory_consumption () const
 

Private Attributes

std::vector< boolblock_mask
 

Friends

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

Detailed Description

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

The semantics of this class are the same as the related ComponentMask class, i.e., a default constructed mask represents all possible blocks. See there for more information about these semantics.

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

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

Note that by wrapping the velocity elements into a single FESystem object we make sure that the overall element has only 2 blocks. The result is a block mask that, in both 2d and 3d, would have values [false, true]. (Compare this to the corresponding component mask discussed in the ComponentMask documentation.) Similarly, using

BlockMask velocity_mask = stokes_fe.block_mask (velocities);

would result in a mask [true, false] in both 2d and 3d.

Definition at line 71 of file block_mask.h.

Constructor & Destructor Documentation

◆ BlockMask() [1/4]

BlockMask::BlockMask ( )
default

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

◆ BlockMask() [2/4]

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

Deprecated constructor allowing implicit conversion.

◆ BlockMask() [3/4]

BlockMask::BlockMask ( const std::vector< bool > &  block_mask)
explicit

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

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

◆ BlockMask() [4/4]

BlockMask::BlockMask ( const unsigned int  n_blocks,
const bool  initializer 
)

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

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

Member Function Documentation

◆ size()

unsigned int BlockMask::size ( ) const

If this block 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 BlockMask::operator[] ( const unsigned int  block_index) const

Return whether a particular block is selected by this mask. If this mask represents the case of an object that selects all blocks (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
block_indexThe index for which the function should return whether the block is selected. If this object represents a mask in which all blocks are always selected then any index is allowed here. Otherwise, the given index needs to be between zero and the number of blocks that this mask represents.

◆ represents_n_blocks()

bool BlockMask::represents_n_blocks ( const unsigned int  n) const

Return whether this block mask represents a mask with exactly n blocks. 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 blocks and will always say that a block is selected.

◆ n_selected_blocks()

unsigned int BlockMask::n_selected_blocks ( const unsigned int  overall_number_of_blocks = numbers::invalid_unsigned_int) const

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

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

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_block()

unsigned int BlockMask::first_selected_block ( const unsigned int  overall_number_of_blocks = numbers::invalid_unsigned_int) const

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

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

◆ represents_the_all_selected_mask()

bool BlockMask::represents_the_all_selected_mask ( ) const

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

◆ operator|()

BlockMask BlockMask::operator| ( const BlockMask mask) const

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

◆ operator&()

BlockMask BlockMask::operator& ( const BlockMask mask) const

Return a block 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 BlockMask::operator== ( const BlockMask mask) const

Return whether this object and the argument are identical.

◆ operator!=()

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

Return whether this object and the argument are not identical.

◆ memory_consumption()

std::size_t BlockMask::memory_consumption ( ) const

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

Definition at line 46 of file block_mask.cc.

Friends And Related Symbol Documentation

◆ operator<<

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

Write a block mask to an output stream. If the block mask represents one where all blocks are selected without specifying a particular size of the mask, then it writes the string [all blocks selected] to the stream. Otherwise, it prints the block 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 block_mask.cc.

Member Data Documentation

◆ block_mask

std::vector<bool> BlockMask::block_mask
private

The actual block mask.

Definition at line 222 of file block_mask.h.


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