Reference documentation for deal.II version 8.5.1
|
#include <deal.II/lac/trilinos_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 std::vector< Epetra_Map > &partitioner) 1 | |
BlockVector (const std::vector< IndexSet > &partitioner, const MPI_Comm &communicator=MPI_COMM_WORLD) 1 | |
BlockVector (const MPI::BlockVector &V) 1 | |
BlockVector (const BlockVector &V) 1 | |
BlockVector (const size_type num_blocks) 1 | |
BlockVector (const std::vector< size_type > &N) 1 | |
template<typename InputIterator > | |
BlockVector (const std::vector< size_type > &n, const InputIterator first, const InputIterator end) 1 | |
~BlockVector () | |
BlockVector & | operator= (const value_type s) |
BlockVector & | operator= (const MPI::BlockVector &V) |
BlockVector & | operator= (const BlockVector &V) |
template<typename Number > | |
BlockVector & | operator= (const ::BlockVector< Number > &V) |
void | reinit (const std::vector< Epetra_Map > &partitioning, const bool omit_zeroing_entries=false) |
void | reinit (const std::vector< IndexSet > &partitioning, const MPI_Comm &communicator=MPI_COMM_WORLD, const bool omit_zeroing_entries=false) |
void | reinit (const std::vector< size_type > &N, const bool omit_zeroing_entries=false) |
void | reinit (const MPI::BlockVector &V) |
void | reinit (const BlockVector &V, const bool omit_zeroing_entries=false) |
void | reinit (const size_type 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 ::ExceptionBase & | ExcNonMatchingBlockVectors () |
static ::ExceptionBase & | ExcNonLocalizedMap (int arg1, int arg2) |
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 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.
In contrast to the class MPI::BlockVector, this class is based on a localized version of the vectors, which means that the whole vector is stored on each processor. Note that matrix vector products with this block vector class do only work in case the program is run on only one processor, since the Trilinos matrices are inherently parallel.
This class is deprecated, use TrilinosWrappers::MPI::BlockVector instead.
Definition at line 70 of file trilinos_block_vector.h.
Typedef the base class for simpler access to its own typedefs.
Definition at line 76 of file trilinos_block_vector.h.
Typedef the type of the underlying vector.
Definition at line 81 of file trilinos_block_vector.h.
typedef BaseClass::value_type TrilinosWrappers::BlockVector::value_type |
Import the typedefs from the base class.
Definition at line 86 of file trilinos_block_vector.h.
|
inline |
Default constructor. Generate an empty vector without any blocks.
Definition at line 320 of file trilinos_block_vector.h.
|
inlineexplicit |
Constructor. Generate a block vector with as many blocks as there are entries in Input_Maps. For this non-distributed vector, the parallel partitioning is not used, just the global size of the partitioner.
Definition at line 326 of file trilinos_block_vector.h.
|
inlineexplicit |
Constructor. Generate a block vector with as many blocks as there are entries in Input_Maps. For this non-distributed vector, the parallel partitioning is not used, just the global size of the partitioner.
Definition at line 334 of file trilinos_block_vector.h.
|
inline |
Copy-Constructor. Set all the properties of the non-parallel vector to those of the given parallel vector and import the elements.
Definition at line 389 of file trilinos_block_vector.h.
|
inline |
Copy-Constructor. Set all the properties of the vector to those of the given input vector and copy the elements.
Definition at line 397 of file trilinos_block_vector.h.
|
inlineexplicit |
Creates a block vector consisting of num_blocks
components, but there is no content in the individual components and the user has to fill appropriate data using a reinit of the blocks.
Definition at line 375 of file trilinos_block_vector.h.
|
inlineexplicit |
Constructor. Set the number of blocks to n.size()
and initialize each block with n[i]
zero elements.
References BlockVector.reinit().
Definition at line 343 of file trilinos_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 351 of file trilinos_block_vector.h.
|
inline |
Destructor. Clears memory
Definition at line 383 of file trilinos_block_vector.h.
BlockVector & BlockVector< Number >::operator= | ( | const value_type | s | ) |
Copy operator: fill all components of the vector that are locally stored with the given scalar value.
Definition at line 261 of file trilinos_block_vector.cc.
BlockVector & BlockVector< Number >::operator= | ( | const MPI::BlockVector & | V | ) |
Copy operator for a distributed Trilinos vector to a localized one.
Definition at line 378 of file trilinos_block_vector.cc.
BlockVector & BlockVector< Number >::operator= | ( | const BlockVector & | V | ) |
Copy operator for arguments of the same type.
Definition at line 388 of file trilinos_block_vector.cc.
BlockVector & BlockVector< Number >::operator= | ( | const ::BlockVector< Number > & | V | ) |
Another copy function. This one takes a deal.II block vector and copies it into a TrilinosWrappers block vector. Note that the number of blocks has to be the same in the vector as in the input vector. Use the reinit() command for resizing the BlockVector or for changing the internal structure of the block components.
Since Trilinos only works on doubles, this function is limited to accept only one possible number type in the deal.II vector.
Definition at line 423 of file trilinos_block_vector.h.
void BlockVector< Number >::reinit | ( | const std::vector< Epetra_Map > & | partitioning, |
const bool | omit_zeroing_entries = false |
||
) |
Reinitialize the BlockVector to contain as many blocks as there are Epetra_Maps given in the input argument, according to the global size of the individual components described in the maps. Note that the resulting vector will be stored completely on each process. The Epetra_Map is useful when data exchange with a distributed vector based on the same Epetra_map is intended. In that case, the same communicator is used for data exchange.
If omit_zeroing_entries==false
, the vector is filled with zeros.
Definition at line 270 of file trilinos_block_vector.cc.
void BlockVector< Number >::reinit | ( | const std::vector< IndexSet > & | partitioning, |
const MPI_Comm & | communicator = MPI_COMM_WORLD , |
||
const bool | omit_zeroing_entries = false |
||
) |
Reinitialize the BlockVector to contain as many blocks as there are index sets given in the input argument, according to the global size of the individual components described in the index set, and using a given MPI communicator. The MPI communicator is useful when data exchange with a distributed vector based on the same initialization is intended. In that case, the same communicator is used for data exchange.
If omit_zeroing_entries==false
, the vector is filled with zeros.
Definition at line 293 of file trilinos_block_vector.cc.
void BlockVector< Number >::reinit | ( | const std::vector< size_type > & | N, |
const bool | omit_zeroing_entries = false |
||
) |
Reinitialize the BlockVector to contain as many blocks as there are elements in the first argument, and with the respective sizes. Since no distribution map is given, all vectors are local vectors.
If omit_zeroing_entries==false
, the vector is filled with zeros.
Definition at line 317 of file trilinos_block_vector.cc.
void BlockVector< Number >::reinit | ( | const MPI::BlockVector & | V | ) |
Reinitialize the vector in the same way as the given to a distributed block vector. The elements will be copied in this process.
Definition at line 333 of file trilinos_block_vector.cc.
void BlockVector< Number >::reinit | ( | const BlockVector & | V, |
const bool | omit_zeroing_entries = false |
||
) |
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 362 of file trilinos_block_vector.cc.
void BlockVector< Number >::reinit | ( | const size_type | 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 346 of file trilinos_block_vector.cc.
|
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 411 of file trilinos_block_vector.h.
void BlockVector< Number >::print | ( | std::ostream & | out, |
const unsigned int | precision = 3 , |
||
const bool | scientific = true , |
||
const bool | across = true |
||
) | const |
Print to a stream.
Definition at line 408 of file trilinos_block_vector.cc.
|
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 451 of file trilinos_block_vector.h.