Reference documentation for deal.II version 9.6.0
|
#include <deal.II/lac/trilinos_parallel_block_vector.h>
Public Types | |
using | BaseClass = ::BlockVectorBase<MPI::Vector> |
using | BlockType = BaseClass::BlockType |
using | value_type = BaseClass::value_type |
using | pointer = BaseClass::pointer |
using | const_pointer = BaseClass::const_pointer |
using | reference = BaseClass::reference |
using | const_reference = BaseClass::const_reference |
using | size_type = BaseClass::size_type |
using | iterator = BaseClass::iterator |
using | const_iterator = BaseClass::const_iterator |
using | real_type |
Public Member Functions | |
BlockVector ()=default | |
BlockVector (const std::vector< IndexSet > ¶llel_partitioning, const MPI_Comm communicator=MPI_COMM_WORLD) | |
BlockVector (const std::vector< IndexSet > ¶llel_partitioning, const std::vector< IndexSet > &ghost_values, const MPI_Comm communicator, const bool vector_writable=false) | |
BlockVector (const BlockVector &v) | |
BlockVector (BlockVector &&v) noexcept | |
BlockVector (const size_type num_blocks) | |
~BlockVector () override=default | |
BlockVector & | operator= (const value_type s) |
BlockVector & | operator= (const BlockVector &v) |
BlockVector & | operator= (BlockVector &&v) noexcept |
template<typename Number > | |
BlockVector & | operator= (const ::BlockVector< Number > &v) |
void | reinit (const std::vector< IndexSet > ¶llel_partitioning, const MPI_Comm communicator=MPI_COMM_WORLD, const bool omit_zeroing_entries=false) |
void | reinit (const std::vector< IndexSet > &partitioning, const std::vector< IndexSet > &ghost_values, const MPI_Comm communicator=MPI_COMM_WORLD, const bool vector_writable=false) |
void | reinit (const std::vector< std::shared_ptr< const Utilities::MPI::Partitioner > > &partitioners, const bool make_ghosted=true, const bool vector_writable=false) |
void | reinit (const BlockVector &V, const bool omit_zeroing_entries=false) |
void | reinit (const size_type num_blocks) |
void | import_nonlocal_data_for_fe (const TrilinosWrappers::BlockSparseMatrix &m, const BlockVector &v) |
bool | has_ghost_elements () const |
void | swap (BlockVector &v) noexcept |
void | print (std::ostream &out, const unsigned int precision=3, const bool scientific=true, const bool across=true) const |
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 |
virtual size_type | size () const override |
std::size_t | locally_owned_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 |
virtual void | extract_subvector_to (const ArrayView< const types::global_dof_index > &indices, ArrayView< value_type > &entries) const override |
void | extract_subvector_to (ForwardIterator indices_begin, const ForwardIterator indices_end, OutputIterator values_begin) const |
virtual void | extract_subvector_to (const ArrayView< const types::global_dof_index > &indices, ArrayView< MPI::Vector::value_type > &elements) const=0 |
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 | update_ghost_values () const |
std::size_t | memory_consumption () const |
template<class Archive > | |
void | serialize (Archive &ar, const unsigned int version) |
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 & | ExcIteratorRangeDoesNotMatchVectorSize () |
static ::ExceptionBase & | ExcNonMatchingBlockVectors () |
static ::ExceptionBase & | ExcInUse (int arg1, std::string arg2, std::string arg3) |
static ::ExceptionBase & | ExcNoSubscriber (std::string arg1, std::string arg2) |
Protected Attributes | |
std::vector< MPI::Vector > | components |
BlockIndices | block_indices |
Private Types | |
using | map_value_type = decltype(counter_map)::value_type |
using | map_iterator = decltype(counter_map)::iterator |
Private Member Functions | |
void | check_no_subscribers () const noexcept |
Private Attributes | |
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 Private Attributes | |
static std::mutex | mutex |
Related Symbols | |
(Note that these are not member symbols.) | |
void | swap (BlockVector &u, BlockVector &v) noexcept |
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.
The model of distribution of data is such that each of the blocks is distributed across all MPI processes named in the MPI communicator. I.e. we don't just distribute the whole vector, but each component. In the constructors and reinit() functions, one therefore not only has to specify the sizes of the individual blocks, but also the number of elements of each of these blocks to be stored on the local process.
Definition at line 74 of file trilinos_parallel_block_vector.h.
Typedef the base class for simpler access to its own alias.
Definition at line 80 of file trilinos_parallel_block_vector.h.
Typedef the type of the underlying vector.
Definition at line 85 of file trilinos_parallel_block_vector.h.
Import the alias from the base class.
Definition at line 90 of file trilinos_parallel_block_vector.h.
Definition at line 91 of file trilinos_parallel_block_vector.h.
Definition at line 92 of file trilinos_parallel_block_vector.h.
Definition at line 93 of file trilinos_parallel_block_vector.h.
Definition at line 94 of file trilinos_parallel_block_vector.h.
Definition at line 95 of file trilinos_parallel_block_vector.h.
Definition at line 96 of file trilinos_parallel_block_vector.h.
Definition at line 97 of file trilinos_parallel_block_vector.h.
|
inherited |
Declare a type that has holds real-valued numbers with the same precision as the template argument to this class. If the template argument of this class is a real data type, then real_type equals the template argument. If the template argument is a std::complex type then real_type equals the type underlying the complex numbers.
This alias is used to represent the return type of norms.
Definition at line 480 of file block_vector_base.h.
|
privateinherited |
The data type used in counter_map.
Definition at line 229 of file subscriptor.h.
|
privateinherited |
The iterator type used in counter_map.
Definition at line 234 of file subscriptor.h.
|
default |
Default constructor. Generate an empty vector without any blocks.
|
inlineexplicit |
Constructor. Generate a block vector with as many blocks as there are entries in partitioning
. Each IndexSet together with the MPI communicator contains the layout of the distribution of data among the MPI processes.
Definition at line 340 of file trilinos_parallel_block_vector.h.
|
inline |
Creates a BlockVector with ghost elements. See the respective reinit() method for more details. ghost_values
may contain any elements in parallel_partitioning
, they will be ignored.
Definition at line 349 of file trilinos_parallel_block_vector.h.
|
inline |
Copy-Constructor. Set all the properties of the parallel vector to those of the given argument and copy the elements.
Definition at line 370 of file trilinos_parallel_block_vector.h.
|
inlinenoexcept |
Move constructor. Creates a new vector by stealing the internal data of the vector v
.
Definition at line 384 of file trilinos_parallel_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 363 of file trilinos_parallel_block_vector.h.
|
overridedefault |
Destructor. Clears memory
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 30 of file trilinos_block_vector.cc.
BlockVector & BlockVector< Number >::operator= | ( | const BlockVector & | v | ) |
Copy operator for arguments of the same type. Resize the present vector if necessary to the correct number of blocks, then copy the individual blocks from v
using the copy-assignment operator of the class that represents the individual blocks.
Copying the vectors that make up individual blocks can have complex semantics in parallel vector classes. See the information provided by the class used to represent the individual blocks.
Definition at line 39 of file trilinos_block_vector.cc.
|
noexcept |
Move the given vector. This operator replaces the present vector with v
by efficiently swapping the internal data structures.
Definition at line 61 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 395 of file trilinos_parallel_block_vector.h.
void BlockVector< Number >::reinit | ( | const std::vector< IndexSet > & | parallel_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 parallel distribution of the individual components described in the maps.
If omit_zeroing_entries==false
, the vector is filled with zeros.
Definition at line 70 of file trilinos_block_vector.cc.
void BlockVector< Number >::reinit | ( | const std::vector< IndexSet > & | partitioning, |
const std::vector< IndexSet > & | ghost_values, | ||
const MPI_Comm | communicator = MPI_COMM_WORLD, | ||
const bool | vector_writable = false ) |
Reinit functionality. This function destroys the old vector content and generates a new one based on the input partitioning. In addition to just specifying one index set as in all the other methods above, this method allows to supply an additional set of ghost entries. There are two different versions of a vector that can be created. If the flag vector_writable
is set to false
, the vector only allows read access to the joint set of parallel_partitioning
and ghost_entries
. The effect of the reinit method is then equivalent to calling the other reinit method with an index set containing both the locally owned entries and the ghost entries.
If the flag vector_writable
is set to true, this creates an alternative storage scheme for ghost elements that allows multiple threads to write into the vector (for the other reinit methods, only one thread is allowed to write into the ghost entries at a time).
Definition at line 91 of file trilinos_block_vector.cc.
void BlockVector< Number >::reinit | ( | const std::vector< std::shared_ptr< const Utilities::MPI::Partitioner > > & | partitioners, |
const bool | make_ghosted = true, | ||
const bool | vector_writable = false ) |
Initialize each block given to each parallel partitioning described in partitioners
.
You can decide whether your vector will contain ghost elements with make_ghosted
.
The parameter vector_writable
only has effect on ghosted vectors and is ignored for non-ghosted vectors.
Definition at line 116 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 139 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 156 of file trilinos_block_vector.cc.
void BlockVector< Number >::import_nonlocal_data_for_fe | ( | const TrilinosWrappers::BlockSparseMatrix & | m, |
const BlockVector & | v ) |
This reinit function is meant to be used for parallel calculations where some non-local data has to be used. The typical situation where one needs this function is the call of the FEValues<dim>::get_function_values function (or of some derivatives) in parallel. Since it is usually faster to retrieve the data in advance, this function can be called before the assembly forks out to the different processors. What this function does is the following: It takes the information in the columns of the given matrix and looks which data couples between the different processors. That data is then queried from the input vector. Note that you should not write to the resulting vector any more, since the some data can be stored several times on different processors, leading to unpredictable results. In particular, such a vector cannot be used for matrix-vector products as for example done during the solution of linear systems.
Definition at line 168 of file trilinos_block_vector.cc.
|
inline |
Return if this Vector contains ghost elements.
Definition at line 417 of file trilinos_parallel_block_vector.h.
|
inlinenoexcept |
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.
Definition at line 430 of file trilinos_parallel_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 189 of file trilinos_block_vector.cc.
|
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 |
Call the compress() function on all the subblocks of the vector.
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.
|
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.
|
overridevirtualinherited |
Return dimension of the vector. This is the sum of the dimensions of all components.
Implements ReadVector< MPI::Vector::value_type >.
|
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 index set that describes which elements of this vector are owned by the current processor. Note that this index set does not include elements this vector may store locally as ghost elements but that are in fact owned by another 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
For block vectors, this function returns the union of the locally owned elements of the individual blocks, shifted by their respective index offsets.
|
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.
|
overridevirtualinherited |
|
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
.
|
pure virtualinherited |
Extract a subset of the vector specified by indices
into the output array elements
.
|
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 |
Return the square of the \(l_2\)-norm.
|
inherited |
Return the mean value of the elements of this vector.
|
inherited |
Return the \(l_1\)-norm of the vector, i.e. the sum of the absolute values.
|
inherited |
Return the \(l_2\)-norm of the vector, i.e. the square root of the sum of the squares of the elements.
|
inherited |
Return the maximum absolute value of the elements of this vector, which is the \(l_\infty\)-norm of a vector.
|
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 whether the vector contains only elements with value zero. This function is mainly for internal consistency check and should seldom be used when not in debug mode since it uses quite some time.
|
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 |
Update the ghost values by calling update_ghost_values
for each block.
|
inherited |
Determine an estimate for the memory consumption (in bytes) of this object.
|
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 135 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 155 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 300 of file subscriptor.h.
|
inlineinherited |
List the subscribers to the input stream
.
Definition at line 317 of file subscriptor.h.
|
inherited |
List the subscribers to deallog
.
Definition at line 203 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 309 of file subscriptor.h.
|
privatenoexceptinherited |
Check that there are no objects subscribing to this object. If this check passes then it is safe to destroy the current object. It this check fails then this function will either abort or print an error message to deallog (by using the AssertNothrow mechanism), but will not throw an exception.
Definition at line 52 of file subscriptor.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 447 of file trilinos_parallel_block_vector.h.
|
protectedinherited |
Pointer to the array of components.
Definition at line 963 of file block_vector_base.h.
|
protectedinherited |
Object managing the transformation between global indices and indices within the different blocks.
Definition at line 969 of file block_vector_base.h.
|
mutableprivateinherited |
Store the number of objects which subscribed to this object. Initially, this number is zero, and upon destruction it shall be zero again (i.e. all objects which subscribed should have unsubscribed again).
The creator (and owner) of an object is counted in the map below if HE manages to supply identification.
We use the mutable
keyword in order to allow subscription to constant objects also.
This counter may be read from and written to concurrently in multithreaded code: hence we use the std::atomic
class template.
Definition at line 218 of file subscriptor.h.
|
mutableprivateinherited |
In this map, we count subscriptions for each different identification string supplied to subscribe().
Definition at line 224 of file subscriptor.h.
|
mutableprivateinherited |
In this vector, we store pointers to the validity bool in the SmartPointer objects that subscribe to this class.
Definition at line 240 of file subscriptor.h.
|
mutableprivateinherited |
Pointer to the typeinfo object of this object, from which we can later deduce the class name. Since this information on the derived class is neither available in the destructor, nor in the constructor, we obtain it in between and store it here.
Definition at line 248 of file subscriptor.h.
|
staticprivateinherited |
A mutex used to ensure data consistency when accessing the mutable
members of this class. This lock is used in the subscribe() and unsubscribe() functions, as well as in list_subscribers()
.
Definition at line 271 of file subscriptor.h.