|
void | collect_sizes () |
|
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 | locally_owned_size () 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 |
|
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 | is_non_negative () const |
|
BlockVectorBase & | operator+= (const BlockVectorBase &V) |
|
BlockVectorBase & | operator-= (const BlockVectorBase &V) |
|
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) |
|
virtual void | compress (VectorOperation::values) |
|
|
| BlockVector (const size_type num_blocks=0, const size_type block_size=0) |
|
| BlockVector (const BlockVector< Number > &V) |
|
template<typename OtherNumber > |
| BlockVector (const BlockVector< OtherNumber > &v) |
|
| BlockVector (const std::vector< size_type > &block_sizes) |
|
| BlockVector (const std::vector< IndexSet > &local_ranges, const std::vector< IndexSet > &ghost_indices, const MPI_Comm &communicator) |
|
| BlockVector (const std::vector< IndexSet > &local_ranges, const MPI_Comm &communicator) |
|
virtual | ~BlockVector () override=default |
|
virtual BlockVector & | operator= (const value_type s) override |
|
BlockVector & | operator= (const BlockVector &V) |
|
template<class Number2 > |
BlockVector & | operator= (const BlockVector< Number2 > &V) |
|
BlockVector & | operator= (const Vector< Number > &V) |
|
BlockVector< Number > & | operator= (const PETScWrappers::MPI::BlockVector &petsc_vec) |
|
BlockVector< Number > & | operator= (const TrilinosWrappers::MPI::BlockVector &trilinos_vec) |
|
void | reinit (const size_type num_blocks, const size_type block_size=0, const bool omit_zeroing_entries=false) |
|
void | reinit (const std::vector< size_type > &block_sizes, const bool omit_zeroing_entries=false) |
|
template<typename Number2 > |
void | reinit (const BlockVector< Number2 > &V, const bool omit_zeroing_entries=false) |
|
virtual void | compress (::VectorOperation::values operation) override |
|
void | update_ghost_values () const |
|
void | zero_out_ghosts () const |
|
void | zero_out_ghost_values () const |
|
bool | has_ghost_elements () const |
|
void | set_ghost_state (const bool ghosted) const |
|
template<typename Number2 > |
void | copy_locally_owned_data_from (const BlockVector< Number2 > &src) |
|
template<typename OtherNumber > |
void | add (const std::vector< size_type > &indices, const ::Vector< OtherNumber > &values) |
|
void | sadd (const Number s, const BlockVector< Number > &V) |
|
virtual bool | all_zero () const override |
|
virtual Number | mean_value () const override |
|
real_type | lp_norm (const real_type p) const |
|
void | swap (BlockVector< Number > &v) |
|
|
virtual void | reinit (const VectorSpaceVector< Number > &V, const bool omit_zeroing_entries=false) override |
|
virtual BlockVector< Number > & | operator*= (const Number factor) override |
|
virtual BlockVector< Number > & | operator/= (const Number factor) override |
|
virtual BlockVector< Number > & | operator+= (const VectorSpaceVector< Number > &V) override |
|
virtual BlockVector< Number > & | operator-= (const VectorSpaceVector< Number > &V) override |
|
virtual void | import (const LinearAlgebra::ReadWriteVector< Number > &V, VectorOperation::values operation, std::shared_ptr< const Utilities::MPI::CommunicationPatternBase > communication_pattern={}) override |
|
virtual Number | operator* (const VectorSpaceVector< Number > &V) const override |
|
template<typename FullMatrixType > |
void | multivector_inner_product (FullMatrixType &matrix, const BlockVector< Number > &V, const bool symmetric=false) const |
|
template<typename FullMatrixType > |
Number | multivector_inner_product_with_metric (const FullMatrixType &matrix, const BlockVector< Number > &V, const bool symmetric=false) const |
|
template<typename FullMatrixType > |
void | mmult (BlockVector< Number > &V, const FullMatrixType &matrix, const Number s=Number(0.), const Number b=Number(1.)) const |
|
virtual void | add (const Number a) override |
|
virtual void | add (const Number a, const VectorSpaceVector< Number > &V) override |
|
virtual void | add (const Number a, const VectorSpaceVector< Number > &V, const Number b, const VectorSpaceVector< Number > &W) override |
|
virtual void | add (const std::vector< size_type > &indices, const std::vector< Number > &values) |
|
virtual void | sadd (const Number s, const Number a, const VectorSpaceVector< Number > &V) override |
|
virtual void | scale (const VectorSpaceVector< Number > &scaling_factors) override |
|
virtual void | equ (const Number a, const VectorSpaceVector< Number > &V) override |
|
virtual real_type | l1_norm () const override |
|
virtual real_type | l2_norm () const override |
|
real_type | norm_sqr () const |
|
virtual real_type | linfty_norm () const override |
|
virtual Number | add_and_dot (const Number a, const VectorSpaceVector< Number > &V, const VectorSpaceVector< Number > &W) override |
|
virtual size_type | size () const override |
|
virtual ::IndexSet | locally_owned_elements () const override |
|
virtual void | print (std::ostream &out, const unsigned int precision=3, const bool scientific=true, const bool across=true) const override |
|
virtual std::size_t | memory_consumption () const override |
|
static ::ExceptionBase & | ExcVectorTypeNotCompatible () |
|
static ::ExceptionBase & | ExcIteratorRangeDoesNotMatchVectorSize () |
|
template<typename Number>
class LinearAlgebra::distributed::BlockVector< Number >
An implementation of block vectors based on distributed deal.II vectors. 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.
- Note
- Instantiations for this template are provided for
<float> and <double>
; others can be generated in application programs (see the section on Template instantiations in the manual).
- See also
- Block (linear algebra)
Definition at line 84 of file la_parallel_block_vector.h.
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 ::Vector< OtherNumber > &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}\).