Reference documentation for deal.II version 9.2.0
\(\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\}}\)
Public Types | Public Member Functions | Static Public Member Functions | Related Functions | List of all members

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

Inheritance diagram for TrilinosWrappers::MPI::BlockVector:
[legend]

Public Types

using BaseClass = ::BlockVectorBase< MPI::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
 
- Public Types inherited from BlockVectorBase< MPI::Vector >
using BlockType = MPI::Vector
 
using value_type = typename BlockType::value_type
 
using pointer = value_type *
 
using const_pointer = const value_type *
 
using iterator = ::internal::BlockVectorIterators::Iterator< BlockVectorBase, false >
 
using const_iterator = ::internal::BlockVectorIterators::Iterator< BlockVectorBase, true >
 
using reference = typename BlockType::reference
 
using const_reference = typename BlockType::const_reference
 
using size_type = types::global_dof_index
 
using real_type = typename BlockType::real_type
 

Public Member Functions

 BlockVector ()=default
 
 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_values, const MPI_Comm &communicator, const bool vector_writable=false)
 
 BlockVector (const BlockVector &v)
 
 BlockVector (BlockVector &&v) noexcept
 
 BlockVector (const size_type num_blocks)
 
 ~BlockVector () override=default
 
BlockVectoroperator= (const value_type s)
 
BlockVectoroperator= (const BlockVector &v)
 
BlockVectoroperator= (BlockVector &&v) noexcept
 
template<typename Number >
BlockVectoroperator= (const ::BlockVector< Number > &v)
 
void reinit (const std::vector< IndexSet > &parallel_partitioning, const MPI_Comm &communicator=MPI_COMM_WORLD, const bool omit_zeroing_entries=false)
 
void reinit (const std::vector< IndexSet > &partitioning, const std::vector< IndexSet > &ghost_values, const MPI_Comm &communicator=MPI_COMM_WORLD, const bool vector_writable=false)
 
void reinit (const BlockVector &V, const bool omit_zeroing_entries=false)
 
void reinit (const size_type num_blocks)
 
void import_nonlocal_data_for_fe (const TrilinosWrappers::BlockSparseMatrix &m, const BlockVector &v)
 
bool has_ghost_elements () 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< MPI::Vector >
 BlockVectorBase ()=default
 
 BlockVectorBase (const BlockVectorBase &)=default
 
 BlockVectorBase (BlockVectorBase &&) noexcept=default
 
void collect_sizes ()
 
void compress (::VectorOperation::values operation)
 
BlockTypeblock (const unsigned int i)
 
const BlockTypeblock (const unsigned int i) const
 
const BlockIndicesget_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
 
BlockVectorBaseoperator= (const value_type s)
 
BlockVectorBaseoperator= (const BlockVectorBase &V)
 
BlockVectorBaseoperator= (BlockVectorBase &&)=default
 
BlockVectorBaseoperator= (const BlockVectorBase< VectorType2 > &V)
 
BlockVectorBaseoperator= (const MPI::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
 
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
 
- Public Member Functions inherited from Subscriptor
 Subscriptor ()
 
 Subscriptor (const Subscriptor &)
 
 Subscriptor (Subscriptor &&) noexcept
 
virtual ~Subscriptor ()
 
Subscriptoroperator= (const Subscriptor &)
 
Subscriptoroperator= (Subscriptor &&) noexcept
 
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
 
template<class Archive >
void serialize (Archive &ar, const unsigned int version)
 

Static Public Member Functions

static ::ExceptionBaseExcIteratorRangeDoesNotMatchVectorSize ()
 
static ::ExceptionBaseExcNonMatchingBlockVectors ()
 
- Static Public Member Functions inherited from Subscriptor
static ::ExceptionBaseExcInUse (int arg1, std::string arg2, std::string arg3)
 
static ::ExceptionBaseExcNoSubscriber (std::string arg1, std::string arg2)
 

Related Functions

(Note that these are not member functions.)

void swap (BlockVector &u, BlockVector &v)
 

Additional Inherited Members

- Protected Attributes inherited from BlockVectorBase< MPI::Vector >
std::vector< MPI::Vector > components
 
BlockIndices block_indices
 

Detailed Description

An implementation of block vectors based on the vector class implemented in TrilinosWrappers. 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)

Author
Martin Kronbichler, Wolfgang Bangerth, 2008, 2009

Definition at line 75 of file trilinos_parallel_block_vector.h.


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