Reference documentation for deal.II version 9.1.1
|
#include <deal.II/lac/block_vector.h>
Public Types | |
using | BaseClass = BlockVectorBase< Vector< Number > > |
using | BlockType = typename BaseClass::BlockType |
using | value_type = typename BaseClass::value_type |
Public Types inherited from BlockVectorBase< Vector< Number > > | |
using | BlockType = Vector< Number > |
using | real_type = typename BlockType::real_type |
Public Member Functions | |
BlockVector (const unsigned int n_blocks=0, const size_type block_size=0) | |
BlockVector (const BlockVector< Number > &V) | |
BlockVector (BlockVector< Number > &&) noexcept=default | |
template<typename OtherNumber > | |
BlockVector (const BlockVector< OtherNumber > &v) | |
BlockVector (const TrilinosWrappers::MPI::BlockVector &v) | |
BlockVector (const std::vector< size_type > &block_sizes) | |
BlockVector (const BlockIndices &block_indices) | |
template<typename InputIterator > | |
BlockVector (const std::vector< size_type > &block_sizes, const InputIterator first, const InputIterator end) | |
~BlockVector () override=default | |
void | compress (::VectorOperation::values operation=::VectorOperation::unknown) |
bool | has_ghost_elements () const |
BlockVector & | operator= (const value_type s) |
BlockVector< Number > & | operator= (const BlockVector< Number > &v) |
BlockVector< Number > & | operator= (BlockVector< Number > &&)=default |
template<class Number2 > | |
BlockVector< Number > & | operator= (const BlockVector< Number2 > &V) |
BlockVector< Number > & | operator= (const Vector< Number > &V) |
BlockVector< Number > & | operator= (const TrilinosWrappers::MPI::BlockVector &V) |
void | reinit (const unsigned int n_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) |
void | reinit (const BlockIndices &block_indices, const bool omit_zeroing_entries=false) |
template<typename Number2 > | |
void | reinit (const BlockVector< Number2 > &V, const bool omit_zeroing_entries=false) |
template<class BlockVector2 > | |
void | scale (const BlockVector2 &v) |
void | swap (BlockVector< Number > &v) |
void | print (std::ostream &out, const unsigned int precision=3, const bool scientific=true, const bool across=true) const |
void | block_write (std::ostream &out) const |
void | block_read (std::istream &in) |
Public Member Functions inherited from BlockVectorBase< Vector< Number > > | |
BlockVectorBase ()=default | |
BlockVectorBase (const BlockVectorBase &)=default | |
BlockVectorBase (BlockVectorBase &&) noexcept=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< Number > &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 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 &&) noexcept | |
virtual | ~Subscriptor () |
Subscriptor & | operator= (const Subscriptor &) |
Subscriptor & | operator= (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 ::ExceptionBase & | ExcIteratorRangeDoesNotMatchVectorSize () |
Static Public Member Functions inherited from Subscriptor | |
static ::ExceptionBase & | ExcInUse (int arg1, std::string arg2, std::string arg3) |
static ::ExceptionBase & | ExcNoSubscriber (std::string arg1, std::string arg2) |
Related Functions | |
(Note that these are not member functions.) | |
template<typename Number > | |
void | swap (BlockVector< Number > &u, BlockVector< Number > &v) |
template<typename Number > | |
void | swap (LinearAlgebra::distributed::BlockVector< Number > &u, LinearAlgebra::distributed::BlockVector< Number > &v) |
Additional Inherited Members | |
Protected Attributes inherited from BlockVectorBase< Vector< Number > > | |
std::vector< Vector< Number > > | components |
BlockIndices | block_indices |
An implementation of block vectors based on 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.
<float> and <double>
; others can be generated in application programs (see the section on Template instantiations in the manual).Definition at line 41 of file block_linear_operator.h.
using BlockVector< Number >::BaseClass = BlockVectorBase<Vector<Number> > |
Typedef the base class for simpler access to its own alias.
Definition at line 74 of file block_vector.h.
using BlockVector< Number >::BlockType = typename BaseClass::BlockType |
Typedef the type of the underlying vector.
Definition at line 79 of file block_vector.h.
using BlockVector< Number >::value_type = typename BaseClass::value_type |
Import the alias from the base class.
Definition at line 84 of file block_vector.h.
|
explicit |
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 n_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.
BlockVector< Number >::BlockVector | ( | const BlockVector< Number > & | V | ) |
Copy Constructor. Dimension set to that of v
, all components are copied from v
.
|
defaultnoexcept |
Move constructor. Creates a new vector by stealing the internal data of the given argument vector.
|
explicit |
Copy constructor taking a BlockVector of another data type. This will fail if there is no conversion path from OtherNumber
to Number
. Note that you may lose accuracy when copying to a BlockVector with data elements with less accuracy.
Older versions of gcc did not honor the explicit
keyword on template constructors. In such cases, it is easy to accidentally write code that can be very inefficient, since the compiler starts performing hidden conversions. To avoid this, this function is disabled if we have detected a broken compiler during configuration.
BlockVector< Number >::BlockVector | ( | const TrilinosWrappers::MPI::BlockVector< Number > & | v | ) |
A copy constructor taking a (parallel) Trilinos block vector and copying it into the deal.II own format.
BlockVector< Number >::BlockVector | ( | const std::vector< size_type > & | block_sizes | ) |
Constructor. Set the number of blocks to block_sizes.size()
and initialize each block with block_sizes[i]
zero elements.
BlockVector< Number >::BlockVector | ( | const BlockIndices & | block_indices | ) |
Constructor. Initialize vector to the structure found in the BlockIndices argument.
BlockVector< Number >::BlockVector | ( | const std::vector< size_type > & | block_sizes, |
const InputIterator | first, | ||
const InputIterator | end | ||
) |
Constructor. Set the number of blocks to block_sizes.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.
|
overridedefault |
Destructor. Clears memory
void BlockVector< Number >::compress | ( | ::VectorOperation::values | operation = ::VectorOperation::unknown | ) |
Call the compress() function on all the subblocks.
This functionality only needs to be called if using MPI based vectors and exists in other objects for compatibility.
See Compressing distributed objects for more information.
bool BlockVector< Number >::has_ghost_elements | ( | ) | const |
Returns false
as this is a serial block vector.
This functionality only needs to be called if using MPI based vectors and exists in other objects for compatibility.
BlockVector& BlockVector< Number >::operator= | ( | const value_type | s | ) |
Copy operator: fill all components of the vector with the given scalar value.
BlockVector<Number>& BlockVector< Number >::operator= | ( | const BlockVector< Number > & | v | ) |
Copy operator for arguments of the same type. Resize the present vector if necessary.
|
default |
Move the given vector. This operator replaces the present vector with the contents of the given argument vector.
BlockVector<Number>& BlockVector< Number >::operator= | ( | const BlockVector< Number2 > & | V | ) |
Copy operator for template arguments of different types. Resize the present vector if necessary.
BlockVector<Number>& BlockVector< Number >::operator= | ( | const Vector< Number > & | V | ) |
Copy a regular vector into a block vector.
BlockVector<Number>& BlockVector< Number >::operator= | ( | const TrilinosWrappers::MPI::BlockVector< Number > & | V | ) |
A copy constructor from a Trilinos block vector to a deal.II block vector.
void BlockVector< Number >::reinit | ( | const unsigned int | n_blocks, |
const size_type | block_size = 0 , |
||
const bool | omit_zeroing_entries = false |
||
) |
Reinitialize the BlockVector to contain n_blocks
blocks of size block_size
each.
If the second argument is left at its default value, then the block vector allocates the specified number of blocks but leaves them at zero size. You then need to later reinitialize the individual blocks, and call collect_sizes() to update the block system's knowledge of its individual block's sizes.
If omit_zeroing_entries==false
, the vector is filled with zeros.
void BlockVector< Number >::reinit | ( | const std::vector< size_type > & | block_sizes, |
const bool | omit_zeroing_entries = false |
||
) |
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.
void BlockVector< Number >::reinit | ( | const BlockIndices & | block_indices, |
const bool | omit_zeroing_entries = false |
||
) |
Reinitialize the BlockVector to reflect the structure found in BlockIndices.
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.
void BlockVector< Number >::reinit | ( | const BlockVector< Number2 > & | 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() of one of the blocks, then subsequent actions of this object may yield unpredictable results since they may be routed to the wrong block.
void BlockVector< Number >::scale | ( | const BlockVector2 & | v | ) |
Multiply each element of this vector by the corresponding element of v
.
void BlockVector< Number >::swap | ( | BlockVector< Number > & | v | ) |
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.
This function is analogous to 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.
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.
void BlockVector< Number >::block_write | ( | std::ostream & | out | ) | const |
Write the vector en bloc to a stream. This is done in a binary mode, so the output is neither readable by humans nor (probably) by other computers using a different operating system or number format.
void BlockVector< Number >::block_read | ( | std::istream & | in | ) |
Read a vector en block from a file. This is done using the inverse operations to the above function, so it is reasonably fast because the bitstream is not interpreted.
The vector is resized if necessary.
A primitive form of error checking is performed which will recognize the bluntest attempts to interpret some data as a vector stored bitwise to a file, but not more.
|
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 488 of file 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 738 of file la_parallel_block_vector.h.