|
| BlockVector ()=default |
|
| 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 > ¶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 () override=default |
|
BlockVector & | operator= (const value_type s) |
|
BlockVector & | operator= (const BlockVector &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 > ¶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 |
|
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 |
|
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 |
|
void | extract_subvector_to (ForwardIterator indices_begin, const ForwardIterator indices_end, OutputIterator values_begin) const |
|
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 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 | update_ghost_values () const |
|
std::size_t | memory_consumption () const |
|
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.
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
IndexSet complete_index_set(const IndexSet::size_type N)
For block vectors, this function returns the union of the locally owned elements of the individual blocks, shifted by their respective index offsets.
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
.
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
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}\).