![]() |
Reference documentation for deal.II version Git 87eb8ae23d 2021-01-27 12:59:39 +0100
|
#include <deal.II/lac/la_parallel_block_vector.h>
Public Types | |
using | BaseClass = BlockVectorBase< Vector< Number > > |
using | BlockType = typename BaseClass::BlockType |
using | value_type = typename BaseClass::value_type |
using | real_type = typename BaseClass::real_type |
using | pointer = typename BaseClass::pointer |
using | const_pointer = typename BaseClass::const_pointer |
using | reference = typename BaseClass::reference |
using | const_reference = typename BaseClass::const_reference |
using | size_type = typename BaseClass::size_type |
using | iterator = typename BaseClass::iterator |
using | const_iterator = typename BaseClass::const_iterator |
Public Member Functions | |
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 |
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 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) |
template<class Archive > | |
void | serialize (Archive &ar, const unsigned int version) |
virtual void | compress (VectorOperation::values) |
1: Basic operations | |
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 > &N, 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 |
bool | has_ghost_elements () const |
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) |
2: Implementation of VectorSpaceVector | |
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 CommunicationPatternBase > communication_pattern=std::shared_ptr< const CommunicationPatternBase >()) 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 |
Subscriptor functionality | |
Classes derived from Subscriptor provide a facility to subscribe to this object. This is mostly used by the SmartPointer class. | |
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 |
Static Public Member Functions | |
static ::ExceptionBase & | ExcVectorTypeNotCompatible () |
static ::ExceptionBase & | ExcIteratorRangeDoesNotMatchVectorSize () |
static ::ExceptionBase & | ExcInUse (int arg1, std::string arg2, std::string arg3) |
static ::ExceptionBase & | ExcNoSubscriber (std::string arg1, std::string arg2) |
Static Public Attributes | |
static constexpr unsigned int | communication_block_size = 20 |
Protected Attributes | |
std::vector< Vector< Number > > | components |
BlockIndices | block_indices |
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.
<float> and <double>
; others can be generated in application programs (see the section on Template instantiations in the manual).Definition at line 84 of file la_parallel_block_vector.h.
using LinearAlgebra::distributed::BlockVector< Number >::BaseClass = BlockVectorBase<Vector<Number> > |
Typedef the base class for simpler access to its own alias.
Definition at line 105 of file la_parallel_block_vector.h.
using LinearAlgebra::distributed::BlockVector< Number >::BlockType = typename BaseClass::BlockType |
Typedef the type of the underlying vector.
Definition at line 110 of file la_parallel_block_vector.h.
using LinearAlgebra::distributed::BlockVector< Number >::value_type = typename BaseClass::value_type |
Import the alias from the base class.
Definition at line 115 of file la_parallel_block_vector.h.
using LinearAlgebra::distributed::BlockVector< Number >::real_type = typename BaseClass::real_type |
Definition at line 116 of file la_parallel_block_vector.h.
using LinearAlgebra::distributed::BlockVector< Number >::pointer = typename BaseClass::pointer |
Definition at line 117 of file la_parallel_block_vector.h.
using LinearAlgebra::distributed::BlockVector< Number >::const_pointer = typename BaseClass::const_pointer |
Definition at line 118 of file la_parallel_block_vector.h.
using LinearAlgebra::distributed::BlockVector< Number >::reference = typename BaseClass::reference |
Definition at line 119 of file la_parallel_block_vector.h.
using LinearAlgebra::distributed::BlockVector< Number >::const_reference = typename BaseClass::const_reference |
Definition at line 120 of file la_parallel_block_vector.h.
using LinearAlgebra::distributed::BlockVector< Number >::size_type = typename BaseClass::size_type |
Definition at line 121 of file la_parallel_block_vector.h.
using LinearAlgebra::distributed::BlockVector< Number >::iterator = typename BaseClass::iterator |
Definition at line 122 of file la_parallel_block_vector.h.
using LinearAlgebra::distributed::BlockVector< Number >::const_iterator = typename BaseClass::const_iterator |
Definition at line 123 of file la_parallel_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 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.
LinearAlgebra::distributed::BlockVector< Number >::BlockVector | ( | const BlockVector< Number > & | V | ) |
Copy-Constructor. Dimension set to that of V, all components are copied from V
|
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.
LinearAlgebra::distributed::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.
LinearAlgebra::distributed::BlockVector< Number >::BlockVector | ( | const std::vector< IndexSet > & | local_ranges, |
const std::vector< IndexSet > & | ghost_indices, | ||
const MPI_Comm & | communicator | ||
) |
Construct a block vector with an IndexSet for the local range and ghost entries for each block.
LinearAlgebra::distributed::BlockVector< Number >::BlockVector | ( | const std::vector< IndexSet > & | local_ranges, |
const MPI_Comm & | communicator | ||
) |
Same as above but the ghost indices are assumed to be empty.
|
overridevirtualdefault |
Destructor.
|
overridevirtual |
Copy operator: fill all components of the vector with the given scalar value.
BlockVector& LinearAlgebra::distributed::BlockVector< Number >::operator= | ( | const BlockVector< Number > & | V | ) |
Copy operator for arguments of the same type. Resize the present vector if necessary.
BlockVector& LinearAlgebra::distributed::BlockVector< Number >::operator= | ( | const BlockVector< Number2 > & | V | ) |
Copy operator for template arguments of different types. Resize the present vector if necessary.
BlockVector& LinearAlgebra::distributed::BlockVector< Number >::operator= | ( | const Vector< Number > & | V | ) |
Copy a regular vector into a block vector.
BlockVector<Number>& LinearAlgebra::distributed::BlockVector< Number >::operator= | ( | const PETScWrappers::MPI::BlockVector< Number > & | petsc_vec | ) |
Copy the content of a PETSc vector into the calling vector. This function assumes that the vectors layouts have already been initialized to match.
This operator is only available if deal.II was configured with PETSc.
BlockVector<Number>& LinearAlgebra::distributed::BlockVector< Number >::operator= | ( | const TrilinosWrappers::MPI::BlockVector< Number > & | trilinos_vec | ) |
Copy the content of a Trilinos vector into the calling vector. This function assumes that the vectors layouts have already been initialized to match.
This operator is only available if deal.II was configured with Trilinos.
void LinearAlgebra::distributed::BlockVector< Number >::reinit | ( | const size_type | num_blocks, |
const size_type | block_size = 0 , |
||
const bool | omit_zeroing_entries = false |
||
) |
Reinitialize the BlockVector to contain num_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 LinearAlgebra::distributed::BlockVector< Number >::reinit | ( | const std::vector< size_type > & | N, |
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 LinearAlgebra::distributed::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.
|
overridevirtual |
This function copies the data that has accumulated in the data buffer for ghost indices to the owning processor. For the meaning of the argument operation
, see the entry on Compressing distributed vectors and matrices in the glossary.
There are two variants for this function. If called with argument VectorOperation::add
adds all the data accumulated in ghost elements to the respective elements on the owning processor and clears the ghost array afterwards. If called with argument VectorOperation::insert
, a set operation is performed. Since setting elements in a vector with ghost elements is ambiguous (as one can set both the element on the ghost site as well as the owning site), this operation makes the assumption that all data is set correctly on the owning processor. Upon call of compress(VectorOperation::insert), all ghost entries are therefore simply zeroed out (using zero_ghost_values()). In debug mode, a check is performed that makes sure that the data set is actually consistent between processors, i.e., whenever a non-zero ghost element is found, it is compared to the value on the owning processor and an exception is thrown if these elements do not agree.
void LinearAlgebra::distributed::BlockVector< Number >::update_ghost_values | ( | ) | const |
Fills the data field for ghost indices with the values stored in the respective positions of the owning processor. This function is needed before reading from ghosts. The function is const
even though ghost data is changed. This is needed to allow functions with a const
vector to perform the data exchange without creating temporaries.
void LinearAlgebra::distributed::BlockVector< Number >::zero_out_ghosts | ( | ) | const |
This method zeros the entries on ghost dofs, but does not touch locally owned DoFs.
After calling this method, read access to ghost elements of the vector is forbidden and an exception is thrown. Only write access to ghost elements is allowed in this state.
bool LinearAlgebra::distributed::BlockVector< Number >::has_ghost_elements | ( | ) | const |
Return if this Vector contains ghost elements.
void LinearAlgebra::distributed::BlockVector< Number >::add | ( | const std::vector< size_type > & | indices, |
const ::Vector< OtherNumber > & | values | ||
) |
This is a collective add operation that adds a whole set of values stored in values
to the vector components specified by indices
.
void LinearAlgebra::distributed::BlockVector< Number >::sadd | ( | const Number | s, |
const BlockVector< Number > & | V | ||
) |
Scaling and simple vector addition, i.e. *this = s*(*this)+V
.
|
overridevirtual |
Return whether the vector contains only elements with value zero. This function is mainly for internal consistency checks and should seldom be used when not in debug mode since it uses quite some time.
Implements LinearAlgebra::VectorSpaceVector< Number >.
|
overridevirtual |
Compute the mean value of all the entries in the vector.
Implements LinearAlgebra::VectorSpaceVector< Number >.
real_type LinearAlgebra::distributed::BlockVector< Number >::lp_norm | ( | const real_type | p | ) | const |
\(l_p\)-norm of the vector. The pth root of the sum of the pth powers of the absolute values of the elements.
void LinearAlgebra::distributed::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.
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 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.
|
overridevirtual |
Change the dimension to that of the vector V. The elements of V are not copied.
Implements LinearAlgebra::VectorSpaceVector< Number >.
|
overridevirtual |
Multiply the entire vector by a fixed factor.
Implements LinearAlgebra::VectorSpaceVector< Number >.
|
overridevirtual |
Divide the entire vector by a fixed factor.
Implements LinearAlgebra::VectorSpaceVector< Number >.
|
overridevirtual |
Add the vector V
to the present one.
Implements LinearAlgebra::VectorSpaceVector< Number >.
|
overridevirtual |
Subtract the vector V
from the present one.
Implements LinearAlgebra::VectorSpaceVector< Number >.
|
overridevirtual |
Import all the elements present in the vector's IndexSet from the input vector V
. VectorOperation::values operation
is used to decide if the elements in V
should be added to the current vector or replace the current elements. The last parameter can be used if the same communication pattern is used multiple times. This can be used to improve performance.
Implements LinearAlgebra::VectorSpaceVector< Number >.
|
overridevirtual |
Return the scalar product of two vectors.
Implements LinearAlgebra::VectorSpaceVector< Number >.
void LinearAlgebra::distributed::BlockVector< Number >::multivector_inner_product | ( | FullMatrixType & | matrix, |
const BlockVector< Number > & | V, | ||
const bool | symmetric = false |
||
) | const |
Calculate the scalar product between each block of this vector and V
and store the result in a full matrix matrix
. This function computes the result by forming \(A_{ij}=U_i \cdot V_j\) where \(U_i\) and \(V_j\) indicate the \(i\)th block (not element!) of \(U\) and the \(j\)th block of \(V\), respectively. If symmetric
is true
, it is assumed that inner product results in a square symmetric matrix and almost half of the scalar products can be avoided.
Obviously, this function can only be used if all blocks of both vectors are of the same size.
Number LinearAlgebra::distributed::BlockVector< Number >::multivector_inner_product_with_metric | ( | const FullMatrixType & | matrix, |
const BlockVector< Number > & | V, | ||
const bool | symmetric = false |
||
) | const |
Calculate the scalar product between each block of this vector and V
using a metric tensor matrix
. This function computes the result of \( \sum_{ij} A^{ij} U_i \cdot V_j\) where \(U_i\) and \(V_j\) indicate the \(i\)th block (not element) of \(U\) and the \(j\)th block of \(V\), respectively. If symmetric
is true
, it is assumed that \(U_i \cdot V_j\) and \(A^{ij}\) are symmetric matrices and almost half of the scalar products can be avoided.
Obviously, this function can only be used if all blocks of both vectors are of the same size.
void LinearAlgebra::distributed::BlockVector< Number >::mmult | ( | BlockVector< Number > & | V, |
const FullMatrixType & | matrix, | ||
const Number | s = Number(0.) , |
||
const Number | b = Number(1.) |
||
) | const |
Set each block of this vector as follows: \(V^i = s V^i + b \sum_{j} U_j A^{ji}\) where \(V^i\) and \(U_j\) indicate the \(i\)th block (not element) of \(V\) and the \(j\)th block of \(U\), respectively.
Obviously, this function can only be used if all blocks of both vectors are of the same size.
|
overridevirtual |
Add a
to all components. Note that a
is a scalar not a vector.
Implements LinearAlgebra::VectorSpaceVector< Number >.
|
overridevirtual |
Simple addition of a multiple of a vector, i.e. *this += a*V
.
Implements LinearAlgebra::VectorSpaceVector< Number >.
|
overridevirtual |
Multiple addition of scaled vectors, i.e. *this += a*V+b*W
.
Implements LinearAlgebra::VectorSpaceVector< Number >.
|
virtual |
A collective add operation: This function adds a whole set of values stored in values
to the vector components specified by indices
.
|
overridevirtual |
Scaling and simple addition of a multiple of a vector, i.e. *this = s*(*this)+a*V
.
Implements LinearAlgebra::VectorSpaceVector< Number >.
|
overridevirtual |
Scale each element of this vector by the corresponding element in the argument. This function is mostly meant to simulate multiplication (and immediate re-assignment) by a diagonal scaling matrix.
Implements LinearAlgebra::VectorSpaceVector< Number >.
|
overridevirtual |
Assignment *this = a*V
.
Implements LinearAlgebra::VectorSpaceVector< Number >.
|
overridevirtual |
Return the l1 norm of the vector (i.e., the sum of the absolute values of all entries among all processors).
Implements LinearAlgebra::VectorSpaceVector< Number >.
|
overridevirtual |
Return the \(l_2\) norm of the vector (i.e., the square root of the sum of the square of all entries among all processors).
Implements LinearAlgebra::VectorSpaceVector< Number >.
real_type LinearAlgebra::distributed::BlockVector< Number >::norm_sqr | ( | ) | const |
Return the square of the \(l_2\) norm of the vector.
|
overridevirtual |
Return the maximum norm of the vector (i.e., the maximum absolute value among all entries and among all processors).
Implements LinearAlgebra::VectorSpaceVector< Number >.
|
overridevirtual |
Perform 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
The reason this function exists is that this operation involves less memory transfer than calling the two functions separately. 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}\).
Implements LinearAlgebra::VectorSpaceVector< Number >.
|
overridevirtual |
Return the global size of the vector, equal to the sum of the number of locally owned indices among all processors.
Implements LinearAlgebra::VectorSpaceVector< Number >.
|
overridevirtual |
Return an index set that describes which elements of this vector are owned by the current processor. As a consequence, the index sets returned on different processors if this is a distributed vector will form disjoint sets that add up to the complete index set. Obviously, if a vector is created on only one processor, then the result would satisfy
Implements LinearAlgebra::VectorSpaceVector< Number >.
|
overridevirtual |
Print the vector to the output stream out
.
Implements LinearAlgebra::VectorSpaceVector< Number >.
|
overridevirtual |
Return the memory consumption of this class in bytes.
Implements LinearAlgebra::VectorSpaceVector< Number >.
|
inherited |
Update internal structures after resizing vectors. Whenever you reinited a block of a block vector, the internal data structures are corrupted. Therefore, you should call this function after all blocks got their new size.
|
inherited |
Access to a single block.
|
inherited |
Read-only access to a single block.
|
inherited |
Return a reference on the object that describes the mapping between block and global indices. The use of this function is highly deprecated and it should vanish in one of the next versions
|
inherited |
Number of blocks.
|
inherited |
Return an iterator pointing to the first element.
|
inherited |
Return an iterator pointing to the first element of a constant block vector.
|
inherited |
Return an iterator pointing to the element past the end.
|
inherited |
Return an iterator pointing to the element past the end of a constant block vector.
|
inherited |
Access components, returns U(i).
|
inherited |
Access components, returns U(i) as a writeable reference.
|
inherited |
Access components, returns U(i).
Exactly the same as operator().
|
inherited |
Access components, returns U(i) as a writeable reference.
Exactly the same as operator().
|
inherited |
Instead of getting individual elements of a vector via operator(), this function allows getting a whole set of elements at once. The indices of the elements to be read are stated in the first argument, the corresponding values are returned in the second.
If the current vector is called v
, then this function is the equivalent to the code
indices
and values
arrays must be identical.
|
inherited |
Instead of getting individual elements of a vector via operator(), this function allows getting a whole set of elements at once. In contrast to the previous function, this function obtains the indices of the elements by dereferencing all elements of the iterator range provided by the first two arguments, and puts the vector values into memory locations obtained by dereferencing a range of iterators starting at the location pointed to by the third argument.
If the current vector is called v
, then this function is the equivalent to the code
values_begin
as there are iterators between indices_begin
and indices_end
.
|
inherited |
Check for equality of two block vector types. This operation is only allowed if the two vectors already have the same block structure.
|
inherited |
\(U = U * V\): scalar product.
|
inherited |
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
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}\).
|
inherited |
Return true if the given global index is in the local range of this processor. Asks the corresponding block.
|
inherited |
Return true
if the vector has no negative entries, i.e. all entries are zero or positive. This function is used, for example, to check whether refinement indicators are really all positive (or zero).
|
inherited |
Addition operator. Fast equivalent to U.add(1, V)
.
|
inherited |
Subtraction operator. Fast equivalent to U.add(-1, V)
.
|
inherited |
A collective add operation: This function adds a whole set of values stored in values
to the vector components specified by indices
.
|
inherited |
This is a second collective add operation. As a difference, this function takes a deal.II vector of values.
|
inherited |
Take an address where n_elements
are stored contiguously and add them into the vector. Handles all cases which are not covered by the other two add()
functions above.
|
inherited |
\(U(0-DIM)+=s\). Addition of s
to all components. Note that s
is a scalar and not a vector.
|
inherited |
U+=a*V. Simple addition of a scaled vector.
|
inherited |
U+=a*V+b*W. Multiple addition of scaled vectors.
|
inherited |
U=s*U+V. Scaling and simple vector addition.
|
inherited |
U=s*U+a*V. Scaling and simple addition.
|
inherited |
U=s*U+a*V+b*W. Scaling and multiple addition.
|
inherited |
U=s*U+a*V+b*W+c*X. Scaling and multiple addition.
|
inherited |
Scale each element of the vector by a constant value.
|
inherited |
Scale each element of the vector by the inverse of the given value.
|
inherited |
Multiply each element of this vector by the corresponding element of v
.
|
inherited |
U=a*V. Assignment.
|
inherited |
Subscribes a user of the object by storing the pointer validity
. The subscriber may be identified by text supplied as identifier
.
Definition at line 136 of file subscriptor.cc.
|
inherited |
Unsubscribes a user from the object.
identifier
and the validity
pointer must be the same as the one supplied to subscribe(). Definition at line 156 of file subscriptor.cc.
|
inlineinherited |
Return the present number of subscriptions to this object. This allows to use this class for reference counted lifetime determination where the last one to unsubscribe also deletes the object.
Definition at line 301 of file subscriptor.h.
|
inlineinherited |
List the subscribers to the input stream
.
Definition at line 318 of file subscriptor.h.
|
inherited |
List the subscribers to deallog
.
Definition at line 204 of file subscriptor.cc.
|
inlineinherited |
Read or write the data of this object to or from a stream for the purpose of serialization using the BOOST serialization library.
This function does not actually serialize any of the member variables of this class. The reason is that what this class stores is only who subscribes to this object, but who does so at the time of storing the contents of this object does not necessarily have anything to do with who subscribes to the object when it is restored. Consequently, we do not want to overwrite the subscribers at the time of restoring, and then there is no reason to write the subscribers out in the first place.
Definition at line 310 of file subscriptor.h.
|
inlinevirtualinherited |
This function does nothing and only exists for backward compatibility.
Definition at line 226 of file vector_space_vector.h.
|
static |
The chunks size to split communication in update_ghost_values() and compress() calls.
Most common MPI implementations will get slow when too many messages/requests are outstanding. Even when messages are small, say 1 kB only, we should collect enough data with communication_block_size
to cover typical infiniband latencies which are around a few microseconds. Sending 20 kB at a throughput of 5 GB/s takes 4 microseconds, so we should arrive at the bandwidth dominated regime then which is good enough.
Definition at line 100 of file la_parallel_block_vector.h.
|
protectedinherited |
Pointer to the array of components.
Definition at line 941 of file block_vector_base.h.
|
protectedinherited |
Object managing the transformation between global indices and indices within the different blocks.
Definition at line 947 of file block_vector_base.h.