Reference documentation for deal.II version GIT relicensing-224-gc660c0d696 2024-03-28 18:40: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 Types | Public Member Functions | Static Public Member Functions | Protected Attributes | Private Types | Private Member Functions | Private Attributes | Static Private Attributes | Related Symbols | List of all members
PETScWrappers::MPI::BlockVector Class Referenceabstract

#include <deal.II/lac/petsc_block_vector.h>

Inheritance diagram for PETScWrappers::MPI::BlockVector:
Inheritance graph
[legend]

Public Types

using BaseClass = BlockVectorBase< Vector >
 
using BlockType = BaseClass::BlockType
 
using value_type = BaseClass::value_type
 
using pointer = BaseClass::pointer
 
using const_pointer = BaseClass::const_pointer
 
using reference = BaseClass::reference
 
using const_reference = BaseClass::const_reference
 
using size_type = BaseClass::size_type
 
using iterator = BaseClass::iterator
 
using const_iterator = BaseClass::const_iterator
 
using real_type = typename BlockType::real_type
 

Public Member Functions

 BlockVector ()
 
 BlockVector (const unsigned int n_blocks, const MPI_Comm communicator, const size_type block_size, const size_type locally_owned_size)
 
 BlockVector (const BlockVector &V)
 
 BlockVector (const std::vector< size_type > &block_sizes, const MPI_Comm communicator, const std::vector< size_type > &local_elements)
 
 BlockVector (const std::vector< IndexSet > &parallel_partitioning, const MPI_Comm communicator=MPI_COMM_WORLD)
 
 BlockVector (const std::vector< IndexSet > &parallel_partitioning, const std::vector< IndexSet > &ghost_indices, const MPI_Comm communicator)
 
 BlockVector (Vec v)
 
template<size_t num_blocks>
 BlockVector (const std::array< Vec, num_blocks > &)
 
 ~BlockVector () override
 
BlockVectoroperator= (const value_type s)
 
BlockVectoroperator= (const BlockVector &V)
 
void reinit (Vec v)
 
void reinit (const unsigned int n_blocks, const MPI_Comm communicator, const size_type block_size, const size_type locally_owned_size, const bool omit_zeroing_entries=false)
 
void reinit (const std::vector< size_type > &block_sizes, const MPI_Comm communicator, const std::vector< size_type > &locally_owned_sizes, const bool omit_zeroing_entries=false)
 
void reinit (const BlockVector &V, const bool omit_zeroing_entries=false)
 
void reinit (const std::vector< IndexSet > &parallel_partitioning, const MPI_Comm communicator)
 
void reinit (const std::vector< IndexSet > &parallel_partitioning, const std::vector< IndexSet > &ghost_entries, const MPI_Comm communicator)
 
void reinit (const std::vector< std::shared_ptr< const Utilities::MPI::Partitioner > > &partitioners, const bool make_ghosted=true)
 
void collect_sizes ()
 
void compress (VectorOperation::values operation)
 
void reinit (const unsigned int num_blocks)
 
bool has_ghost_elements () const
 
MPI_Comm get_mpi_communicator () const
 
 operator const Vec & () const
 
Vec & petsc_vector ()
 
void swap (BlockVector &v)
 
void print (std::ostream &out, const unsigned int precision=3, const bool scientific=true, const bool across=true) const
 
BlockTypeblock (const unsigned int i)
 
const BlockTypeblock (const unsigned int i) const
 
const BlockIndicesget_block_indices () const
 
unsigned int n_blocks () const
 
virtual size_type size () const override
 
std::size_t locally_owned_size () const
 
IndexSet locally_owned_elements () const
 
iterator begin ()
 
const_iterator begin () const
 
iterator end ()
 
const_iterator end () const
 
value_type operator() (const size_type i) const
 
reference operator() (const size_type i)
 
value_type operator[] (const size_type i) const
 
reference operator[] (const size_type i)
 
void extract_subvector_to (const std::vector< size_type > &indices, std::vector< OtherNumber > &values) const
 
virtual void extract_subvector_to (const ArrayView< const types::global_dof_index > &indices, ArrayView< value_type > &entries) const override
 
void extract_subvector_to (ForwardIterator indices_begin, const ForwardIterator indices_end, OutputIterator values_begin) const
 
virtual void extract_subvector_to (const ArrayView< const types::global_dof_index > &indices, ArrayView< VectorType::value_type > &elements) const=0
 
bool operator== (const BlockVectorBase< VectorType2 > &v) const
 
value_type operator* (const BlockVectorBase &V) const
 
real_type norm_sqr () const
 
value_type mean_value () const
 
real_type l1_norm () const
 
real_type l2_norm () const
 
real_type linfty_norm () const
 
value_type add_and_dot (const value_type a, const BlockVectorBase &V, const BlockVectorBase &W)
 
bool in_local_range (const size_type global_index) const
 
bool all_zero () const
 
bool is_non_negative () const
 
BlockVectorBaseoperator+= (const BlockVectorBase &V)
 
BlockVectorBaseoperator-= (const BlockVectorBase &V)
 
void add (const std::vector< size_type > &indices, const std::vector< Number > &values)
 
void add (const std::vector< size_type > &indices, const Vector< Number > &values)
 
void add (const size_type n_elements, const size_type *indices, const Number *values)
 
void add (const value_type s)
 
void add (const value_type a, const BlockVectorBase &V)
 
void add (const value_type a, const BlockVectorBase &V, const value_type b, const BlockVectorBase &W)
 
void sadd (const value_type s, const BlockVectorBase &V)
 
void sadd (const value_type s, const value_type a, const BlockVectorBase &V)
 
void sadd (const value_type s, const value_type a, const BlockVectorBase &V, const value_type b, const BlockVectorBase &W)
 
void sadd (const value_type s, const value_type a, const BlockVectorBase &V, const value_type b, const BlockVectorBase &W, const value_type c, const BlockVectorBase &X)
 
BlockVectorBaseoperator*= (const value_type factor)
 
BlockVectorBaseoperator/= (const value_type factor)
 
void scale (const BlockVector2 &v)
 
void equ (const value_type a, const BlockVector2 &V)
 
void update_ghost_values () const
 
std::size_t memory_consumption () const
 
template<class Archive >
void serialize (Archive &ar, const unsigned int version)
 
Subscriptor functionality

Classes derived from Subscriptor provide a facility to subscribe to this object. This is mostly used by the SmartPointer class.

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
 

Static Public Member Functions

static ::ExceptionBaseExcIteratorRangeDoesNotMatchVectorSize ()
 
static ::ExceptionBaseExcNonMatchingBlockVectors ()
 
static ::ExceptionBaseExcInUse (int arg1, std::string arg2, std::string arg3)
 
static ::ExceptionBaseExcNoSubscriber (std::string arg1, std::string arg2)
 

Protected Attributes

std::vector< Vectorcomponents
 
BlockIndices block_indices
 

Private Types

using map_value_type = decltype(counter_map)::value_type
 
using map_iterator = decltype(counter_map)::iterator
 

Private Member Functions

void setup_nest_vec ()
 
void check_no_subscribers () const noexcept
 

Private Attributes

Vec petsc_nest_vector
 
std::atomic< unsigned intcounter
 
std::map< std::string, unsigned intcounter_map
 
std::vector< std::atomic< bool > * > validity_pointers
 
const std::type_info * object_info
 

Static Private Attributes

static std::mutex mutex
 

Related Symbols

(Note that these are not member symbols.)

void swap (BlockVector &u, BlockVector &v)
 

Detailed Description

An implementation of block vectors based on the parallel vector class implemented in PETScWrappers. While the base class provides for most of the interface, this class handles the actual allocation of vectors and provides functions that are specific to the underlying vector type.

The model of distribution of data is such that each of the blocks is distributed across all MPI processes named in the MPI communicator. I.e. we don't just distribute the whole vector, but each component. In the constructors and reinit() functions, one therefore not only has to specify the sizes of the individual blocks, but also the number of elements of each of these blocks to be stored on the local process.

@ Block (linear algebra)

Definition at line 60 of file petsc_block_vector.h.

Member Typedef Documentation

◆ BaseClass

Typedef the base class for simpler access to its own alias.

Definition at line 66 of file petsc_block_vector.h.

◆ BlockType

Typedef the type of the underlying vector.

Definition at line 71 of file petsc_block_vector.h.

◆ value_type

Import the alias from the base class.

Definition at line 76 of file petsc_block_vector.h.

◆ pointer

Definition at line 77 of file petsc_block_vector.h.

◆ const_pointer

Definition at line 78 of file petsc_block_vector.h.

◆ reference

Definition at line 79 of file petsc_block_vector.h.

◆ const_reference

Definition at line 80 of file petsc_block_vector.h.

◆ size_type

Definition at line 81 of file petsc_block_vector.h.

◆ iterator

Definition at line 82 of file petsc_block_vector.h.

◆ const_iterator

Definition at line 83 of file petsc_block_vector.h.

◆ real_type

using BlockVectorBase< Vector >::real_type = typename BlockType::real_type
inherited

Declare a type that has holds real-valued numbers with the same precision as the template argument to this class. If the template argument of this class is a real data type, then real_type equals the template argument. If the template argument is a std::complex type then real_type equals the type underlying the complex numbers.

This alias is used to represent the return type of norms.

Definition at line 480 of file block_vector_base.h.

◆ map_value_type

using Subscriptor::map_value_type = decltype(counter_map)::value_type
privateinherited

The data type used in counter_map.

Definition at line 229 of file subscriptor.h.

◆ map_iterator

using Subscriptor::map_iterator = decltype(counter_map)::iterator
privateinherited

The iterator type used in counter_map.

Definition at line 234 of file subscriptor.h.

Constructor & Destructor Documentation

◆ BlockVector() [1/8]

BlockVector< Number >::BlockVector ( )
inline

Default constructor. Generate an empty vector without any blocks.

Definition at line 376 of file petsc_block_vector.h.

◆ BlockVector() [2/8]

BlockVector< Number >::BlockVector ( const unsigned int  n_blocks,
const MPI_Comm  communicator,
const size_type  block_size,
const size_type  locally_owned_size 
)
inlineexplicit

Constructor. Generate a block vector with n_blocks blocks, each of which is a parallel vector across communicator with block_size elements of which locally_owned_size elements are stored on the present process.

Definition at line 383 of file petsc_block_vector.h.

◆ BlockVector() [3/8]

BlockVector< Number >::BlockVector ( const BlockVector V)
inline

Copy constructor. Set all the properties of the parallel vector to those of the given argument and copy the elements.

Definition at line 405 of file petsc_block_vector.h.

◆ BlockVector() [4/8]

BlockVector< Number >::BlockVector ( const std::vector< size_type > &  block_sizes,
const MPI_Comm  communicator,
const std::vector< size_type > &  local_elements 
)
inline

Constructor. Set the number of blocks to block_sizes.size() and initialize each block with block_sizes[i] zero elements. The individual blocks are distributed across the given communicator, and each store local_elements[i] elements on the present process.

Definition at line 394 of file petsc_block_vector.h.

◆ BlockVector() [5/8]

BlockVector< Number >::BlockVector ( const std::vector< IndexSet > &  parallel_partitioning,
const MPI_Comm  communicator = MPI_COMM_WORLD 
)
inlineexplicit

Create a BlockVector with parallel_partitioning.size() blocks, each initialized with the given IndexSet.

Definition at line 419 of file petsc_block_vector.h.

◆ BlockVector() [6/8]

BlockVector< Number >::BlockVector ( const std::vector< IndexSet > &  parallel_partitioning,
const std::vector< IndexSet > &  ghost_indices,
const MPI_Comm  communicator 
)
inline

Same as above, but include ghost elements

Definition at line 429 of file petsc_block_vector.h.

◆ BlockVector() [7/8]

BlockVector< Number >::BlockVector ( Vec  v)
inlineexplicit

Create a BlockVector with a PETSc Vec It infers the number of blocks from the Vec if it is of type VECNEST, otherwise the block vector will only have a single block. Internally, we always store a VECNEST vector.

Definition at line 440 of file petsc_block_vector.h.

◆ BlockVector() [8/8]

template<size_t num_blocks>
BlockVector< num_blocks >::BlockVector ( const std::array< Vec, num_blocks > &  arrayV)
inlineexplicit

Create a BlockVector with an array of PETSc vectors.

Definition at line 449 of file petsc_block_vector.h.

◆ ~BlockVector()

BlockVector< Number >::~BlockVector ( )
override

Destructor. Clears memory

Definition at line 37 of file petsc_parallel_block_vector.cc.

Member Function Documentation

◆ operator=() [1/2]

BlockVector & BlockVector< Number >::operator= ( const value_type  s)
inline

Copy operator: fill all components of the vector that are locally stored with the given scalar value.

Definition at line 463 of file petsc_block_vector.h.

◆ operator=() [2/2]

BlockVector & BlockVector< Number >::operator= ( const BlockVector V)
inline

Copy operator for arguments of the same type.

Definition at line 472 of file petsc_block_vector.h.

◆ reinit() [1/8]

void BlockVector< Number >::reinit ( Vec  v)

This method associates the PETSc Vec to the instance of the class. Infers the number of blocks from v if it is of type VECNEST, otherwise the block vector will only have a single block.

Definition at line 58 of file petsc_parallel_block_vector.cc.

◆ reinit() [2/8]

void BlockVector< Number >::reinit ( const unsigned int  n_blocks,
const MPI_Comm  communicator,
const size_type  block_size,
const size_type  locally_owned_size,
const bool  omit_zeroing_entries = false 
)
inline

Reinitialize the BlockVector to contain n_blocks of size block_size, each of which stores locally_owned_size elements locally. The communicator argument denotes which MPI channel each of these blocks shall communicate.

If omit_zeroing_entries==false, the vector is filled with zeros.

Definition at line 494 of file petsc_block_vector.h.

◆ reinit() [3/8]

void BlockVector< Number >::reinit ( const std::vector< size_type > &  block_sizes,
const MPI_Comm  communicator,
const std::vector< size_type > &  locally_owned_sizes,
const bool  omit_zeroing_entries = false 
)
inline

Reinitialize the BlockVector such that it contains block_sizes.size() blocks. Each block is reinitialized to dimension block_sizes[i]. Each of them stores locally_owned_sizes[i] elements on the present process.

If the number of blocks is the same as before this function was called, all vectors remain the same and reinit() is called for each vector.

If omit_zeroing_entries==false, the vector is filled with zeros.

Note that you must call this (or the other reinit() functions) function, rather than calling the reinit() functions of an individual block, to allow the block vector to update its caches of vector sizes. If you call reinit() of one of the blocks, then subsequent actions on this object may yield unpredictable results since they may be routed to the wrong block.

Definition at line 509 of file petsc_block_vector.h.

◆ reinit() [4/8]

void BlockVector< Number >::reinit ( const BlockVector V,
const bool  omit_zeroing_entries = false 
)
inline

Change the dimension to that of the vector V. The same applies as for the other reinit() function.

The elements of V are not copied, i.e. this function is the same as calling reinit (V.size(), omit_zeroing_entries).

Note that you must call this (or the other reinit() functions) function, rather than calling the reinit() functions of an individual block, to allow the block vector to update its caches of vector sizes. If you call reinit() on one of the blocks, then subsequent actions on this object may yield unpredictable results since they may be routed to the wrong block.

Definition at line 528 of file petsc_block_vector.h.

◆ reinit() [5/8]

void BlockVector< Number >::reinit ( const std::vector< IndexSet > &  parallel_partitioning,
const MPI_Comm  communicator 
)
inline

Reinitialize the BlockVector using IndexSets. See the constructor with the same arguments for details.

Definition at line 543 of file petsc_block_vector.h.

◆ reinit() [6/8]

void BlockVector< Number >::reinit ( const std::vector< IndexSet > &  parallel_partitioning,
const std::vector< IndexSet > &  ghost_entries,
const MPI_Comm  communicator 
)
inline

Same as above but include ghost entries.

Definition at line 561 of file petsc_block_vector.h.

◆ reinit() [7/8]

void BlockVector< Number >::reinit ( const std::vector< std::shared_ptr< const Utilities::MPI::Partitioner > > &  partitioners,
const bool  make_ghosted = true 
)
inline

Initialize each block given to each parallel partitioning described in partitioners.

You can decide whether your vector will contain ghost elements with make_ghosted.

Definition at line 584 of file petsc_block_vector.h.

◆ collect_sizes()

void BlockVector< Number >::collect_sizes ( )

This function collects the sizes of the sub-objects and stores them in internal arrays, in order to be able to relay global indices into the vector to indices into the subobjects. You must call this function each time after you have changed the size of the sub-objects.

Definition at line 121 of file petsc_parallel_block_vector.cc.

◆ compress()

void BlockVector< Number >::compress ( VectorOperation::values  operation)

Call the compress() function on all the subblocks of the vector and update the internal state of the nested PETSc vector.

See Compressing distributed objects for more information.

Definition at line 128 of file petsc_parallel_block_vector.cc.

◆ reinit() [8/8]

void BlockVector< Number >::reinit ( const unsigned int  num_blocks)

Change the number of blocks to num_blocks. The individual blocks will get initialized with zero size, so it is assumed that the user resizes the individual blocks by herself in an appropriate way, and calls collect_sizes afterwards.

Definition at line 44 of file petsc_parallel_block_vector.cc.

◆ has_ghost_elements()

bool BlockVector< Number >::has_ghost_elements ( ) const
inline

Return if this vector is a ghosted vector (and thus read-only).

Definition at line 612 of file petsc_block_vector.h.

◆ get_mpi_communicator()

MPI_Comm BlockVector< Number >::get_mpi_communicator ( ) const
inline

Return the underlying MPI communicator.

Definition at line 604 of file petsc_block_vector.h.

◆ operator const Vec &()

BlockVector< Number >::operator const Vec & ( ) const

Conversion operator to gain access to the underlying PETSc type. If you do this, you cut this class off some information it may need, so this conversion operator should only be used if you know what you do. In particular, it should only be used for read-only operations into the vector.

Definition at line 115 of file petsc_parallel_block_vector.cc.

◆ petsc_vector()

Vec & BlockVector< Number >::petsc_vector ( )

Return a reference to the underlying PETSc type. It can be used to modify the underlying data, so use it only when you know what you are doing.

Definition at line 110 of file petsc_parallel_block_vector.cc.

◆ swap()

void BlockVector< Number >::swap ( BlockVector v)
inline

Swap the contents of this vector and the other vector v. One could do this operation with a temporary variable and copying over the data elements, but this function is significantly more efficient since it only swaps the pointers to the data of the two vectors and therefore does not need to allocate temporary storage and move data around.

Limitation: right now this function only works if both vectors have the same number of blocks. If needed, the numbers of blocks should be exchanged, too.

This function is analogous to the swap() function of all C++ standard containers. Also, there is a global function swap(u,v) that simply calls u.swap(v), again in analogy to standard functions.

Definition at line 625 of file petsc_block_vector.h.

◆ print()

void BlockVector< Number >::print ( std::ostream &  out,
const unsigned int  precision = 3,
const bool  scientific = true,
const bool  across = true 
) const
inline

Print to a stream.

Definition at line 636 of file petsc_block_vector.h.

◆ setup_nest_vec()

void BlockVector< Number >::setup_nest_vec ( )
private

Utility to set up the VECNEST object

Definition at line 135 of file petsc_parallel_block_vector.cc.

◆ block() [1/2]

BlockType & BlockVectorBase< Vector >::block ( const unsigned int  i)
inherited

Access to a single block.

◆ block() [2/2]

const BlockType & BlockVectorBase< Vector >::block ( const unsigned int  i) const
inherited

Read-only access to a single block.

◆ get_block_indices()

const BlockIndices & BlockVectorBase< Vector >::get_block_indices ( ) const
inherited

Return a reference on the object that describes the mapping between block and global indices. The use of this function is highly deprecated and it should vanish in one of the next versions

◆ n_blocks()

unsigned int BlockVectorBase< Vector >::n_blocks ( ) const
inherited

Number of blocks.

◆ size()

virtual size_type BlockVectorBase< Vector >::size ( ) const
overridevirtualinherited

Return dimension of the vector. This is the sum of the dimensions of all components.

Implements ReadVector< VectorType::value_type >.

◆ locally_owned_size()

std::size_t BlockVectorBase< Vector >::locally_owned_size ( ) const
inherited

Return local dimension of the vector. This is the sum of the local dimensions (i.e., values stored on the current processor) of all components.

◆ locally_owned_elements()

IndexSet BlockVectorBase< Vector >::locally_owned_elements ( ) const
inherited

Return an index set that describes which elements of this vector are owned by the current processor. Note that this index set does not include elements this vector may store locally as ghost elements but that are in fact owned by another processor. As a consequence, the index sets returned on different processors if this is a distributed vector will form disjoint sets that add up to the complete index set. Obviously, if a vector is created on only one processor, then the result would satisfy

vec.locally_owned_elements() == complete_index_set (vec.size())
IndexSet complete_index_set(const IndexSet::size_type N)
Definition index_set.h:1193

For block vectors, this function returns the union of the locally owned elements of the individual blocks, shifted by their respective index offsets.

◆ begin() [1/2]

iterator BlockVectorBase< Vector >::begin ( )
inherited

Return an iterator pointing to the first element.

◆ begin() [2/2]

const_iterator BlockVectorBase< Vector >::begin ( ) const
inherited

Return an iterator pointing to the first element of a constant block vector.

◆ end() [1/2]

iterator BlockVectorBase< Vector >::end ( )
inherited

Return an iterator pointing to the element past the end.

◆ end() [2/2]

const_iterator BlockVectorBase< Vector >::end ( ) const
inherited

Return an iterator pointing to the element past the end of a constant block vector.

◆ operator()() [1/2]

value_type BlockVectorBase< Vector >::operator() ( const size_type  i) const
inherited

Access components, returns U(i).

◆ operator()() [2/2]

reference BlockVectorBase< Vector >::operator() ( const size_type  i)
inherited

Access components, returns U(i) as a writeable reference.

◆ operator[]() [1/2]

value_type BlockVectorBase< Vector >::operator[] ( const size_type  i) const
inherited

Access components, returns U(i).

Exactly the same as operator().

◆ operator[]() [2/2]

reference BlockVectorBase< Vector >::operator[] ( const size_type  i)
inherited

Access components, returns U(i) as a writeable reference.

Exactly the same as operator().

◆ extract_subvector_to() [1/4]

void BlockVectorBase< Vector >::extract_subvector_to ( const std::vector< size_type > &  indices,
std::vector< OtherNumber > &  values 
) const
inherited

Instead of getting individual elements of a vector via operator(), this function allows getting a whole set of elements at once. The indices of the elements to be read are stated in the first argument, the corresponding values are returned in the second.

If the current vector is called v, then this function is the equivalent to the code

for (unsigned int i=0; i<indices.size(); ++i)
values[i] = v[indices[i]];
Precondition
The sizes of the indices and values arrays must be identical.

◆ extract_subvector_to() [2/4]

virtual void BlockVectorBase< Vector >::extract_subvector_to ( const ArrayView< const types::global_dof_index > &  indices,
ArrayView< value_type > &  entries 
) const
overridevirtualinherited

◆ extract_subvector_to() [3/4]

void BlockVectorBase< Vector >::extract_subvector_to ( ForwardIterator  indices_begin,
const ForwardIterator  indices_end,
OutputIterator  values_begin 
) const
inherited

Instead of getting individual elements of a vector via operator(), this function allows getting a whole set of elements at once. In contrast to the previous function, this function obtains the indices of the elements by dereferencing all elements of the iterator range provided by the first two arguments, and puts the vector values into memory locations obtained by dereferencing a range of iterators starting at the location pointed to by the third argument.

If the current vector is called v, then this function is the equivalent to the code

ForwardIterator indices_p = indices_begin;
OutputIterator values_p = values_begin;
while (indices_p != indices_end)
{
*values_p = v[*indices_p];
++indices_p;
++values_p;
}
Precondition
It must be possible to write into as many memory locations starting at values_begin as there are iterators between indices_begin and indices_end.

◆ extract_subvector_to() [4/4]

virtual void ReadVector< VectorType::value_type >::extract_subvector_to ( const ArrayView< const types::global_dof_index > &  indices,
ArrayView< VectorType::value_type > &  elements 
) const
pure virtualinherited

Extract a subset of the vector specified by indices into the output array elements.

◆ operator==()

bool BlockVectorBase< Vector >::operator== ( const BlockVectorBase< VectorType2 > &  v) const
inherited

Check for equality of two block vector types. This operation is only allowed if the two vectors already have the same block structure.

◆ operator*()

value_type BlockVectorBase< Vector >::operator* ( const BlockVectorBase< Vector > &  V) const
inherited

\(U = U * V\): scalar product.

◆ norm_sqr()

real_type BlockVectorBase< Vector >::norm_sqr ( ) const
inherited

Return the square of the \(l_2\)-norm.

◆ mean_value()

value_type BlockVectorBase< Vector >::mean_value ( ) const
inherited

Return the mean value of the elements of this vector.

◆ l1_norm()

real_type BlockVectorBase< Vector >::l1_norm ( ) const
inherited

Return the \(l_1\)-norm of the vector, i.e. the sum of the absolute values.

◆ l2_norm()

real_type BlockVectorBase< Vector >::l2_norm ( ) const
inherited

Return the \(l_2\)-norm of the vector, i.e. the square root of the sum of the squares of the elements.

◆ linfty_norm()

real_type BlockVectorBase< Vector >::linfty_norm ( ) const
inherited

Return the maximum absolute value of the elements of this vector, which is the \(l_\infty\)-norm of a vector.

◆ add_and_dot()

value_type BlockVectorBase< Vector >::add_and_dot ( const value_type  a,
const BlockVectorBase< Vector > &  V,
const BlockVectorBase< Vector > &  W 
)
inherited

Performs a combined operation of a vector addition and a subsequent inner product, returning the value of the inner product. In other words, the result of this function is the same as if the user called

this->add(a, V);
return_value = *this * W;
void add(const std::vector< size_type > &indices, const std::vector< Number > &values)

The reason this function exists is that this operation involves less memory transfer than calling the two functions separately on deal.II's vector classes (Vector<Number> and LinearAlgebra::distributed::Vector<double>). This method only needs to load three vectors, this, V, W, whereas calling separate methods means to load the calling vector this twice. Since most vector operations are memory transfer limited, this reduces the time by 25% (or 50% if W equals this).

For complex-valued vectors, the scalar product in the second step is implemented as \(\left<v,w\right>=\sum_i v_i \bar{w_i}\).

◆ in_local_range()

bool BlockVectorBase< Vector >::in_local_range ( const size_type  global_index) const
inherited

Return true if the given global index is in the local range of this processor. Asks the corresponding block.

◆ all_zero()

bool BlockVectorBase< Vector >::all_zero ( ) const
inherited

Return whether the vector contains only elements with value zero. This function is mainly for internal consistency check and should seldom be used when not in debug mode since it uses quite some time.

◆ is_non_negative()

bool BlockVectorBase< Vector >::is_non_negative ( ) const
inherited

Return true if the vector has no negative entries, i.e. all entries are zero or positive. This function is used, for example, to check whether refinement indicators are really all positive (or zero).

◆ operator+=()

BlockVectorBase & BlockVectorBase< Vector >::operator+= ( const BlockVectorBase< Vector > &  V)
inherited

Addition operator. Fast equivalent to U.add(1, V).

◆ operator-=()

BlockVectorBase & BlockVectorBase< Vector >::operator-= ( const BlockVectorBase< Vector > &  V)
inherited

Subtraction operator. Fast equivalent to U.add(-1, V).

◆ add() [1/6]

void BlockVectorBase< Vector >::add ( const std::vector< size_type > &  indices,
const std::vector< Number > &  values 
)
inherited

A collective add operation: This function adds a whole set of values stored in values to the vector components specified by indices.

◆ add() [2/6]

void BlockVectorBase< Vector >::add ( const std::vector< size_type > &  indices,
const Vector< Number > &  values 
)
inherited

This is a second collective add operation. As a difference, this function takes a deal.II vector of values.

◆ add() [3/6]

void BlockVectorBase< Vector >::add ( const size_type  n_elements,
const size_type indices,
const Number *  values 
)
inherited

Take an address where n_elements are stored contiguously and add them into the vector. Handles all cases which are not covered by the other two add() functions above.

◆ add() [4/6]

void BlockVectorBase< Vector >::add ( const value_type  s)
inherited

\(U(0-DIM)+=s\). Addition of s to all components. Note that s is a scalar and not a vector.

◆ add() [5/6]

void BlockVectorBase< Vector >::add ( const value_type  a,
const BlockVectorBase< Vector > &  V 
)
inherited

U+=a*V. Simple addition of a scaled vector.

◆ add() [6/6]

void BlockVectorBase< Vector >::add ( const value_type  a,
const BlockVectorBase< Vector > &  V,
const value_type  b,
const BlockVectorBase< Vector > &  W 
)
inherited

U+=a*V+b*W. Multiple addition of scaled vectors.

◆ sadd() [1/4]

void BlockVectorBase< Vector >::sadd ( const value_type  s,
const BlockVectorBase< Vector > &  V 
)
inherited

U=s*U+V. Scaling and simple vector addition.

◆ sadd() [2/4]

void BlockVectorBase< Vector >::sadd ( const value_type  s,
const value_type  a,
const BlockVectorBase< Vector > &  V 
)
inherited

U=s*U+a*V. Scaling and simple addition.

◆ sadd() [3/4]

void BlockVectorBase< Vector >::sadd ( const value_type  s,
const value_type  a,
const BlockVectorBase< Vector > &  V,
const value_type  b,
const BlockVectorBase< Vector > &  W 
)
inherited

U=s*U+a*V+b*W. Scaling and multiple addition.

◆ sadd() [4/4]

void BlockVectorBase< Vector >::sadd ( const value_type  s,
const value_type  a,
const BlockVectorBase< Vector > &  V,
const value_type  b,
const BlockVectorBase< Vector > &  W,
const value_type  c,
const BlockVectorBase< Vector > &  X 
)
inherited

U=s*U+a*V+b*W+c*X. Scaling and multiple addition.

◆ operator*=()

BlockVectorBase & BlockVectorBase< Vector >::operator*= ( const value_type  factor)
inherited

Scale each element of the vector by a constant value.

◆ operator/=()

BlockVectorBase & BlockVectorBase< Vector >::operator/= ( const value_type  factor)
inherited

Scale each element of the vector by the inverse of the given value.

◆ scale()

void BlockVectorBase< Vector >::scale ( const BlockVector2 &  v)
inherited

Multiply each element of this vector by the corresponding element of v.

◆ equ()

void BlockVectorBase< Vector >::equ ( const value_type  a,
const BlockVector2 &  V 
)
inherited

U=a*V. Assignment.

◆ update_ghost_values()

void BlockVectorBase< Vector >::update_ghost_values ( ) const
inherited

Update the ghost values by calling update_ghost_values for each block.

◆ memory_consumption()

std::size_t BlockVectorBase< Vector >::memory_consumption ( ) const
inherited

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

◆ subscribe()

void Subscriptor::subscribe ( std::atomic< bool > *const  validity,
const std::string &  identifier = "" 
) const
inherited

Subscribes a user of the object by storing the pointer validity. The subscriber may be identified by text supplied as identifier.

Definition at line 135 of file subscriptor.cc.

◆ unsubscribe()

void Subscriptor::unsubscribe ( std::atomic< bool > *const  validity,
const std::string &  identifier = "" 
) const
inherited

Unsubscribes a user from the object.

Note
The identifier and the validity pointer must be the same as the one supplied to subscribe().

Definition at line 155 of file subscriptor.cc.

◆ n_subscriptions()

unsigned int Subscriptor::n_subscriptions ( ) const
inlineinherited

Return the present number of subscriptions to this object. This allows to use this class for reference counted lifetime determination where the last one to unsubscribe also deletes the object.

Definition at line 300 of file subscriptor.h.

◆ list_subscribers() [1/2]

template<typename StreamType >
void Subscriptor::list_subscribers ( StreamType &  stream) const
inlineinherited

List the subscribers to the input stream.

Definition at line 317 of file subscriptor.h.

◆ list_subscribers() [2/2]

void Subscriptor::list_subscribers ( ) const
inherited

List the subscribers to deallog.

Definition at line 203 of file subscriptor.cc.

◆ serialize()

template<class Archive >
void Subscriptor::serialize ( Archive &  ar,
const unsigned int  version 
)
inlineinherited

Read or write the data of this object to or from a stream for the purpose of serialization using the BOOST serialization library.

This function does not actually serialize any of the member variables of this class. The reason is that what this class stores is only who subscribes to this object, but who does so at the time of storing the contents of this object does not necessarily have anything to do with who subscribes to the object when it is restored. Consequently, we do not want to overwrite the subscribers at the time of restoring, and then there is no reason to write the subscribers out in the first place.

Definition at line 309 of file subscriptor.h.

◆ check_no_subscribers()

void Subscriptor::check_no_subscribers ( ) const
privatenoexceptinherited

Check that there are no objects subscribing to this object. If this check passes then it is safe to destroy the current object. It this check fails then this function will either abort or print an error message to deallog (by using the AssertNothrow mechanism), but will not throw an exception.

Note
Since this function is just a consistency check it does nothing in release mode.
If this function is called when there is an uncaught exception then, rather than aborting, this function prints an error message to the standard error stream and returns.

Definition at line 52 of file subscriptor.cc.

Friends And Related Symbol Documentation

◆ swap()

void swap ( BlockVector u,
BlockVector v 
)
related

Global function which overloads the default implementation of the C++ standard library which uses a temporary object. The function simply exchanges the data of the two vectors.

Definition at line 661 of file petsc_block_vector.h.

Member Data Documentation

◆ petsc_nest_vector

Vec PETScWrappers::MPI::BlockVector::petsc_nest_vector
private

A PETSc Vec object that describes the entire block vector. Internally, this is done by creating a "nested" vector using PETSc's VECNEST object whose individual blocks are the blocks of this vector.

Definition at line 363 of file petsc_block_vector.h.

◆ components

std::vector<Vector > BlockVectorBase< Vector >::components
protectedinherited

Pointer to the array of components.

Definition at line 963 of file block_vector_base.h.

◆ block_indices

BlockIndices BlockVectorBase< Vector >::block_indices
protectedinherited

Object managing the transformation between global indices and indices within the different blocks.

Definition at line 969 of file block_vector_base.h.

◆ counter

std::atomic<unsigned int> Subscriptor::counter
mutableprivateinherited

Store the number of objects which subscribed to this object. Initially, this number is zero, and upon destruction it shall be zero again (i.e. all objects which subscribed should have unsubscribed again).

The creator (and owner) of an object is counted in the map below if HE manages to supply identification.

We use the mutable keyword in order to allow subscription to constant objects also.

This counter may be read from and written to concurrently in multithreaded code: hence we use the std::atomic class template.

Definition at line 218 of file subscriptor.h.

◆ counter_map

std::map<std::string, unsigned int> Subscriptor::counter_map
mutableprivateinherited

In this map, we count subscriptions for each different identification string supplied to subscribe().

Definition at line 224 of file subscriptor.h.

◆ validity_pointers

std::vector<std::atomic<bool> *> Subscriptor::validity_pointers
mutableprivateinherited

In this vector, we store pointers to the validity bool in the SmartPointer objects that subscribe to this class.

Definition at line 240 of file subscriptor.h.

◆ object_info

const std::type_info* Subscriptor::object_info
mutableprivateinherited

Pointer to the typeinfo object of this object, from which we can later deduce the class name. Since this information on the derived class is neither available in the destructor, nor in the constructor, we obtain it in between and store it here.

Definition at line 248 of file subscriptor.h.

◆ mutex

std::mutex Subscriptor::mutex
staticprivateinherited

A mutex used to ensure data consistency when accessing the mutable members of this class. This lock is used in the subscribe() and unsubscribe() functions, as well as in list_subscribers().

Definition at line 271 of file subscriptor.h.


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