Reference documentation for deal.II version 8.5.1
|
#include <deal.II/lac/petsc_parallel_block_vector.h>
Public Types | |
typedef BlockVectorBase< Vector > | BaseClass |
typedef BaseClass::BlockType | BlockType |
typedef BaseClass::value_type | value_type |
Public Types inherited from BlockVectorBase< Vector > | |
typedef Vector | BlockType |
typedef BlockType::real_type | real_type |
Public Member Functions | |
BlockVector () | |
BlockVector (const unsigned int n_blocks, const MPI_Comm &communicator, const size_type block_size, const size_type local_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 > ¶llel_partitioning, const MPI_Comm &communicator=MPI_COMM_WORLD) | |
BlockVector (const std::vector< IndexSet > ¶llel_partitioning, const std::vector< IndexSet > &ghost_indices, const MPI_Comm &communicator) | |
~BlockVector () | |
BlockVector & | operator= (const value_type s) |
BlockVector & | operator= (const BlockVector &V) |
BlockVector & | operator= (const PETScWrappers::BlockVector &v) |
void | reinit (const unsigned int n_blocks, const MPI_Comm &communicator, const size_type block_size, const size_type local_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 > &local_sizes, const bool omit_zeroing_entries=false) |
void | reinit (const BlockVector &V, const bool omit_zeroing_entries=false) |
void | reinit (const std::vector< IndexSet > ¶llel_partitioning, const MPI_Comm &communicator) |
void | reinit (const std::vector< IndexSet > ¶llel_partitioning, const std::vector< IndexSet > &ghost_entries, const MPI_Comm &communicator) |
void | reinit (const unsigned int num_blocks) |
bool | has_ghost_elements () const |
const MPI_Comm & | get_mpi_communicator () const |
void | swap (BlockVector &v) |
void | print (std::ostream &out, const unsigned int precision=3, const bool scientific=true, const bool across=true) const |
Public Member Functions inherited from BlockVectorBase< Vector > | |
BlockVectorBase () | |
BlockVectorBase (const BlockVectorBase &V)=default | |
BlockVectorBase (BlockVectorBase &&)=default | |
void | collect_sizes () |
void | compress (::VectorOperation::values operation) |
BlockType & | block (const unsigned int i) |
const BlockType & | block (const unsigned int i) const |
const BlockIndices & | get_block_indices () const |
unsigned int | n_blocks () const |
std::size_t | 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 |
void | extract_subvector_to (ForwardIterator indices_begin, const ForwardIterator indices_end, OutputIterator values_begin) const |
BlockVectorBase & | operator= (const value_type s) |
BlockVectorBase & | operator= (const BlockVectorBase &V) |
BlockVectorBase & | operator= (BlockVectorBase &&)=default |
BlockVectorBase & | operator= (const BlockVectorBase< VectorType2 > &V) |
BlockVectorBase & | operator= (const Vector &v) |
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 |
BlockVectorBase & | operator+= (const BlockVectorBase &V) |
BlockVectorBase & | operator-= (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 BlockVectorBase &V) 1 |
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) |
BlockVectorBase & | operator*= (const value_type factor) |
BlockVectorBase & | operator/= (const value_type factor) |
void | scale (const BlockVector2 &v) |
void | equ (const value_type a, const BlockVector2 &V) |
void | equ (const value_type a, const BlockVectorBase &V, const value_type b, const BlockVectorBase &W) |
void | update_ghost_values () const |
std::size_t | memory_consumption () const |
Public Member Functions inherited from Subscriptor | |
Subscriptor () | |
Subscriptor (const Subscriptor &) | |
Subscriptor (Subscriptor &&) | |
virtual | ~Subscriptor () |
Subscriptor & | operator= (const Subscriptor &) |
Subscriptor & | operator= (Subscriptor &&) |
void | subscribe (const char *identifier=0) const |
void | unsubscribe (const char *identifier=0) const |
unsigned int | n_subscriptions () const |
void | list_subscribers () const |
template<class Archive > | |
void | serialize (Archive &ar, const unsigned int version) |
Static Public Member Functions | |
static ::ExceptionBase & | ExcIteratorRangeDoesNotMatchVectorSize () |
static ::ExceptionBase & | ExcNonMatchingBlockVectors () |
Static Public Member Functions inherited from Subscriptor | |
static ::ExceptionBase & | ExcInUse (int arg1, char *arg2, std::string &arg3) |
static ::ExceptionBase & | ExcNoSubscriber (char *arg1, char *arg2) |
Related Functions | |
(Note that these are not member functions.) | |
void | swap (BlockVector &u, BlockVector &v) |
Additional Inherited Members | |
Static Public Attributes inherited from BlockVectorBase< Vector > | |
static const bool | supports_distributed_data |
Protected Attributes inherited from BlockVectorBase< Vector > | |
std::vector< Vector > | components |
BlockIndices | block_indices |
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.
Definition at line 62 of file petsc_parallel_block_vector.h.
Typedef the base class for simpler access to its own typedefs.
Definition at line 68 of file petsc_parallel_block_vector.h.
Typedef the type of the underlying vector.
Definition at line 73 of file petsc_parallel_block_vector.h.
typedef BaseClass::value_type PETScWrappers::MPI::BlockVector::value_type |
Import the typedefs from the base class.
Definition at line 78 of file petsc_parallel_block_vector.h.
|
inline |
Default constructor. Generate an empty vector without any blocks.
Definition at line 305 of file petsc_parallel_block_vector.h.
|
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 local_size
elements are stored on the present process.
Definition at line 311 of file petsc_parallel_block_vector.h.
|
inline |
Copy constructor. Set all the properties of the parallel vector to those of the given argument and copy the elements.
Definition at line 331 of file petsc_parallel_block_vector.h.
|
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 322 of file petsc_parallel_block_vector.h.
|
inlineexplicit |
Create a BlockVector with parallel_partitioning.size() blocks, each initialized with the given IndexSet.
Definition at line 343 of file petsc_parallel_block_vector.h.
|
inline |
Same as above, but include ghost elements
Definition at line 350 of file petsc_parallel_block_vector.h.
|
inline |
Destructor. Clears memory
Definition at line 386 of file petsc_parallel_block_vector.h.
|
inline |
Copy operator: fill all components of the vector that are locally stored with the given scalar value.
Definition at line 359 of file petsc_parallel_block_vector.h.
|
inline |
Copy operator for arguments of the same type.
Definition at line 367 of file petsc_parallel_block_vector.h.
BlockVector & BlockVector< Number >::operator= | ( | const PETScWrappers::BlockVector & | v | ) |
Copy the given sequential (non-distributed) block vector into the present parallel block vector. It is assumed that they have the same size, and this operation does not change the partitioning of the parallel vectors by which its elements are distributed across several MPI processes. What this operation therefore does is to copy that chunk of the given vector v
that corresponds to elements of the target vector that are stored locally, and copies them, for each of the individual blocks of this object. Elements that are not stored locally are not touched.
This being a parallel vector, you must make sure that all processes call this function at the same time. It is not possible to change the local part of a parallel vector on only one process, independent of what other processes do, with this function.
Definition at line 31 of file petsc_parallel_block_vector.cc.
|
inline |
Reinitialize the BlockVector to contain n_blocks
of size block_size
, each of which stores local_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 392 of file petsc_parallel_block_vector.h.
|
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 local_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 408 of file petsc_parallel_block_vector.h.
|
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 425 of file petsc_parallel_block_vector.h.
|
inline |
Reinitialize the BlockVector using IndexSets. See the constructor with the same arguments for details.
Definition at line 438 of file petsc_parallel_block_vector.h.
|
inline |
Same as above but include ghost entries.
Definition at line 455 of file petsc_parallel_block_vector.h.
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 46 of file petsc_parallel_block_vector.cc.
|
inline |
Return if this vector is a ghosted vector (and thus read-only).
Definition at line 482 of file petsc_parallel_block_vector.h.
|
inline |
Return a reference to the MPI communicator object in use with this vector.
Definition at line 475 of file petsc_parallel_block_vector.h.
|
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 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 495 of file petsc_parallel_block_vector.h.
|
inline |
Print to a stream.
Definition at line 506 of file petsc_parallel_block_vector.h.
|
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 532 of file petsc_parallel_block_vector.h.