![]() |
Reference documentation for deal.II version 9.4.1
|
#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 |
std::size_t | locally_owned_size () 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 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) |
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 > &block_sizes, 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 |
void | zero_out_ghost_values () const |
bool | has_ghost_elements () const |
void | set_ghost_state (const bool ghosted) const |
template<typename Number2 > | |
void | copy_locally_owned_data_from (const BlockVector< Number2 > &src) |
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) |
Static Public Attributes | |
static constexpr unsigned int | communication_block_size = 20 |
Protected Attributes | |
std::vector< Vector< Number > > | components |
BlockIndices | block_indices |
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 Utilities::MPI::CommunicationPatternBase > communication_pattern={}) 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 |
static ::ExceptionBase & | ExcVectorTypeNotCompatible () |
static ::ExceptionBase & | ExcIteratorRangeDoesNotMatchVectorSize () |
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 |
template<class Archive > | |
void | serialize (Archive &ar, const unsigned int version) |
static ::ExceptionBase & | ExcInUse (int arg1, std::string arg2, std::string arg3) |
static ::ExceptionBase & | ExcNoSubscriber (std::string arg1, std::string arg2) |
using | map_value_type = decltype(counter_map)::value_type |
using | map_iterator = decltype(counter_map)::iterator |
std::atomic< unsigned int > | counter |
std::map< std::string, unsigned int > | counter_map |
std::vector< std::atomic< bool > * > | validity_pointers |
const std::type_info * | object_info |
static std::mutex | mutex |
void | check_no_subscribers () const noexcept |
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.
|
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 local dimension of the vector. This is the sum of the local dimensions (i.e., values stored on the current processor) of all components.
|
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 |
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.
|
protectedinherited |
Pointer to the array of components.
Definition at line 950 of file block_vector_base.h.
|
protectedinherited |
Object managing the transformation between global indices and indices within the different blocks.
Definition at line 956 of file block_vector_base.h.