Reference documentation for deal.II version 8.5.1
|
#include <deal.II/lac/petsc_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 (const unsigned int num_blocks=0, const size_type block_size=0) | |
BlockVector (const BlockVector &V) | |
BlockVector (const MPI::BlockVector &v) | |
BlockVector (const std::vector< size_type > &n) | |
template<typename InputIterator > | |
BlockVector (const std::vector< size_type > &n, const InputIterator first, const InputIterator end) | |
~BlockVector () | |
BlockVector & | operator= (const value_type s) |
BlockVector & | operator= (const BlockVector &V) |
BlockVector & | operator= (const MPI::BlockVector &v) |
void | reinit (const unsigned int num_blocks, const size_type block_size, const bool omit_zeroing_entries=false) |
void | reinit (const std::vector< size_type > &N, const bool omit_zeroing_entries=false) |
void | reinit (const BlockVector &V, const bool omit_zeroing_entries=false) |
void | reinit (const unsigned int num_blocks) |
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 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 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.
This class is deprecated, use PETScWrappers::MPI::BlockVector.
Definition at line 55 of file petsc_block_vector.h.
Typedef the base class for simpler access to its own typedefs.
Definition at line 61 of file petsc_block_vector.h.
Typedef the type of the underlying vector.
Definition at line 66 of file petsc_block_vector.h.
typedef BaseClass::value_type PETScWrappers::BlockVector::value_type |
Import the typedefs from the base class.
Definition at line 71 of file petsc_block_vector.h.
|
inlineexplicit |
Constructor. There are three ways to use this constructor. First, without any arguments, it generates an object with no blocks. Given one argument, it initializes num_blocks
blocks, but these blocks have size zero. The third variant finally initializes all blocks to the same size block_size
.
Confer the other constructor further down if you intend to use blocks of different sizes.
Definition at line 260 of file petsc_block_vector.h.
|
inline |
Copy-Constructor. Dimension set to that of V, all components are copied from V
Definition at line 276 of file petsc_block_vector.h.
|
inlineexplicit |
Copy-constructor: copy the values from a PETSc wrapper parallel block vector class.
Note that due to the communication model of MPI, all processes have to actually perform this operation, even if they do not use the result. It is not sufficient if only one processor tries to copy the elements from the other processors over to its own process space.
Definition at line 290 of file petsc_block_vector.h.
|
inline |
Constructor. Set the number of blocks to n.size()
and initialize each block with n[i]
zero elements.
Definition at line 269 of file petsc_block_vector.h.
BlockVector< InputIterator >::BlockVector | ( | const std::vector< size_type > & | n, |
const InputIterator | first, | ||
const InputIterator | end | ||
) |
Constructor. Set the number of blocks to n.size()
. Initialize the vector with the elements pointed to by the range of iterators given as second and third argument. Apart from the first argument, this constructor is in complete analogy to the respective constructor of the std::vector
class, but the first argument is needed in order to know how to subdivide the block vector into different blocks.
Definition at line 304 of file petsc_block_vector.h.
|
inline |
Destructor. Clears memory
Definition at line 358 of file petsc_block_vector.h.
|
inline |
Copy operator: fill all components of the vector with the given scalar value.
Definition at line 329 of file petsc_block_vector.h.
|
inline |
Copy operator for arguments of the same type.
Definition at line 339 of file petsc_block_vector.h.
|
inline |
Copy all the elements of the parallel block vector v
into this local vector. Note that due to the communication model of MPI, all processes have to actually perform this operation, even if they do not use the result. It is not sufficient if only one processor tries to copy the elements from the other processors over to its own process space.
Definition at line 349 of file petsc_block_vector.h.
|
inline |
Reinitialize the BlockVector to contain num_blocks
blocks of size block_size
each.
If omit_zeroing_entries==false
, the vector is filled with zeros.
Definition at line 364 of file petsc_block_vector.h.
|
inline |
Reinitialize the BlockVector such that it contains block_sizes.size()
blocks. Each block is reinitialized to dimension block_sizes[i]
.
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() 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 376 of file petsc_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() of one of the blocks, then subsequent actions of this object may yield unpredictable results since they may be routed to the wrong block.
Definition at line 390 of file petsc_block_vector.h.
|
inline |
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 405 of file petsc_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 414 of file petsc_block_vector.h.
|
inline |
Print to a stream.
Definition at line 428 of file petsc_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 455 of file petsc_block_vector.h.