deal.II version GIT relicensing-2165-gc91f007519 2024-11-20 01:40:00+00:00
|
#include <deal.II/lac/block_vector_base.h>
Public Types | |
using | BlockType = VectorType |
using | value_type = typename BlockType::value_type |
using | pointer = value_type * |
using | const_pointer = const value_type * |
using | iterator = ::internal::BlockVectorIterators::Iterator< BlockVectorBase, false > |
using | const_iterator = ::internal::BlockVectorIterators::Iterator< BlockVectorBase, true > |
using | reference = typename BlockType::reference |
using | const_reference = typename BlockType::const_reference |
using | size_type = types::global_dof_index |
using | real_type = typename BlockType::real_type |
Public Member Functions | |
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 |
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) |
template<typename OtherNumber > | |
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, const ArrayView< value_type > &entries) const override |
template<typename ForwardIterator , typename OutputIterator > | |
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 |
template<typename VectorType2 > | |
BlockVectorBase & | operator= (const BlockVectorBase< VectorType2 > &V) |
BlockVectorBase & | operator= (const VectorType &v) |
template<typename VectorType2 > | |
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) |
template<typename Number > | |
void | add (const std::vector< size_type > &indices, const std::vector< Number > &values) |
template<typename Number > | |
void | add (const std::vector< size_type > &indices, const Vector< Number > &values) |
template<typename Number > | |
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) |
template<class BlockVector2 > | |
void | scale (const BlockVector2 &v) |
template<class BlockVector2 > | |
void | equ (const value_type a, const BlockVector2 &V) |
void | update_ghost_values () const |
MPI_Comm | get_mpi_communicator () const |
std::size_t | memory_consumption () const |
virtual void | extract_subvector_to (const ArrayView< const types::global_dof_index > &indices, const ArrayView< VectorType::value_type > &elements) const=0 |
template<class Archive > | |
void | serialize (Archive &ar, const unsigned int version) |
EnableObserverPointer functionality | |
Classes derived from EnableObserverPointer provide a facility to subscribe to this object. This is mostly used by the ObserverPointer 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 & | ExcInUse (int arg1, std::string arg2, std::string arg3) |
static ::ExceptionBase & | ExcNoSubscriber (std::string arg1, std::string arg2) |
Protected Attributes | |
std::vector< VectorType > | 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 |
Friends | |
template<typename N , bool C> | |
class | ::internal::BlockVectorIterators::Iterator |
template<typename > | |
class | BlockVectorBase |
A vector composed of several blocks each representing a vector of its own.
The BlockVector is a collection of vectors of a given type (e.g., deal.II Vector objects, PETScWrappers::MPI::Vector objects, etc.). Each of the vectors inside can have a different size.
The functionality of BlockVector includes everything a Vector can do, plus the access to a single Vector inside the BlockVector by block(i)
. It also has a complete random access iterator, just as the other Vector classes or the standard C++ library template std::vector
. Therefore, all algorithms working on iterators also work with objects of this class.
While this base class implements most of the functionality by dispatching calls to its member functions to the respective functions on each of the individual blocks, this class does not actually allocate some memory or change the size of vectors. For this, the constructors, assignment operators and reinit() functions of derived classes are responsible. This class only handles the common part that is independent of the actual vector type the block vector is built on.
Apart from using this object as a whole, you can use each block separately as a vector, using the block() function. There is a single caveat: if you have changed the size of one or several blocks, you must call the function collect_sizes() of the block vector to update its internal structures.
Definition at line 440 of file block_vector_base.h.
using BlockVectorBase< VectorType >::BlockType = VectorType |
Typedef the type of the underlying vector.
Definition at line 446 of file block_vector_base.h.
using BlockVectorBase< VectorType >::value_type = typename BlockType::value_type |
Definition at line 458 of file block_vector_base.h.
using BlockVectorBase< VectorType >::pointer = value_type * |
Definition at line 459 of file block_vector_base.h.
using BlockVectorBase< VectorType >::const_pointer = const value_type * |
Definition at line 460 of file block_vector_base.h.
using BlockVectorBase< VectorType >::iterator = ::internal::BlockVectorIterators::Iterator<BlockVectorBase, false> |
Definition at line 461 of file block_vector_base.h.
using BlockVectorBase< VectorType >::const_iterator = ::internal::BlockVectorIterators::Iterator<BlockVectorBase, true> |
Definition at line 463 of file block_vector_base.h.
using BlockVectorBase< VectorType >::reference = typename BlockType::reference |
Definition at line 465 of file block_vector_base.h.
using BlockVectorBase< VectorType >::const_reference = typename BlockType::const_reference |
Definition at line 466 of file block_vector_base.h.
using BlockVectorBase< VectorType >::size_type = types::global_dof_index |
Definition at line 467 of file block_vector_base.h.
using BlockVectorBase< VectorType >::real_type = typename BlockType::real_type |
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 478 of file block_vector_base.h.
|
privateinherited |
The data type used in counter_map.
Definition at line 238 of file enable_observer_pointer.h.
|
privateinherited |
The iterator type used in counter_map.
Definition at line 243 of file enable_observer_pointer.h.
|
default |
Default constructor.
|
default |
Copy constructor.
|
defaultnoexcept |
Move constructor. Each block of the argument vector is moved into the current object if the underlying VectorType
is move-constructible, otherwise they are copied.
void BlockVectorBase< VectorType >::collect_sizes | ( | ) |
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.
void BlockVectorBase< VectorType >::compress | ( | VectorOperation::values | operation | ) |
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.
BlockType & BlockVectorBase< VectorType >::block | ( | const unsigned int | i | ) |
Access to a single block.
const BlockType & BlockVectorBase< VectorType >::block | ( | const unsigned int | i | ) | const |
Read-only access to a single block.
const BlockIndices & BlockVectorBase< VectorType >::get_block_indices | ( | ) | const |
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
unsigned int BlockVectorBase< VectorType >::n_blocks | ( | ) | const |
Number of blocks.
|
overridevirtual |
Return dimension of the vector. This is the sum of the dimensions of all components.
Implements ReadVector< VectorType::value_type >.
Reimplemented in LinearAlgebra::distributed::BlockVector< Number >.
std::size_t BlockVectorBase< VectorType >::locally_owned_size | ( | ) | const |
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.
IndexSet BlockVectorBase< VectorType >::locally_owned_elements | ( | ) | const |
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.
iterator BlockVectorBase< VectorType >::begin | ( | ) |
Return an iterator pointing to the first element.
const_iterator BlockVectorBase< VectorType >::begin | ( | ) | const |
Return an iterator pointing to the first element of a constant block vector.
iterator BlockVectorBase< VectorType >::end | ( | ) |
Return an iterator pointing to the element past the end.
const_iterator BlockVectorBase< VectorType >::end | ( | ) | const |
Return an iterator pointing to the element past the end of a constant block vector.
value_type BlockVectorBase< VectorType >::operator() | ( | const size_type | i | ) | const |
Access components, returns U(i).
reference BlockVectorBase< VectorType >::operator() | ( | const size_type | i | ) |
Access components, returns U(i) as a writeable reference.
value_type BlockVectorBase< VectorType >::operator[] | ( | const size_type | i | ) | const |
Access components, returns U(i).
Exactly the same as operator().
reference BlockVectorBase< VectorType >::operator[] | ( | const size_type | i | ) |
Access components, returns U(i) as a writeable reference.
Exactly the same as operator().
void BlockVectorBase< VectorType >::extract_subvector_to | ( | const std::vector< size_type > & | indices, |
std::vector< OtherNumber > & | values | ||
) | const |
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.
|
overridevirtual |
void BlockVectorBase< VectorType >::extract_subvector_to | ( | ForwardIterator | indices_begin, |
const ForwardIterator | indices_end, | ||
OutputIterator | values_begin | ||
) | const |
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
. BlockVectorBase & BlockVectorBase< VectorType >::operator= | ( | const value_type | s | ) |
Copy operator: fill all components of the vector with the given scalar value.
BlockVectorBase & BlockVectorBase< VectorType >::operator= | ( | const BlockVectorBase< VectorType > & | V | ) |
Copy operator for arguments of the same type.
|
default |
Move assignment operator. Move each block of the given argument vector into the current object if VectorType
is move-constructible, otherwise copy them.
BlockVectorBase & BlockVectorBase< VectorType >::operator= | ( | const BlockVectorBase< VectorType2 > & | V | ) |
Copy operator for template arguments of different types.
BlockVectorBase & BlockVectorBase< VectorType >::operator= | ( | const VectorType & | v | ) |
Copy operator from non-block vectors to block vectors.
bool BlockVectorBase< VectorType >::operator== | ( | const BlockVectorBase< VectorType2 > & | v | ) | const |
Check for equality of two block vector types. This operation is only allowed if the two vectors already have the same block structure.
value_type BlockVectorBase< VectorType >::operator* | ( | const BlockVectorBase< VectorType > & | V | ) | const |
\(U = U * V\): scalar product.
real_type BlockVectorBase< VectorType >::norm_sqr | ( | ) | const |
Return the square of the \(l_2\)-norm.
value_type BlockVectorBase< VectorType >::mean_value | ( | ) | const |
Return the mean value of the elements of this vector.
real_type BlockVectorBase< VectorType >::l1_norm | ( | ) | const |
Return the \(l_1\)-norm of the vector, i.e. the sum of the absolute values.
real_type BlockVectorBase< VectorType >::l2_norm | ( | ) | const |
Return the \(l_2\)-norm of the vector, i.e. the square root of the sum of the squares of the elements.
real_type BlockVectorBase< VectorType >::linfty_norm | ( | ) | const |
Return the maximum absolute value of the elements of this vector, which is the \(l_\infty\)-norm of a vector.
value_type BlockVectorBase< VectorType >::add_and_dot | ( | const value_type | a, |
const BlockVectorBase< VectorType > & | V, | ||
const BlockVectorBase< VectorType > & | W | ||
) |
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}\).
bool BlockVectorBase< VectorType >::in_local_range | ( | const size_type | global_index | ) | const |
Return true if the given global index is in the local range of this processor. Asks the corresponding block.
bool BlockVectorBase< VectorType >::all_zero | ( | ) | const |
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.
bool BlockVectorBase< VectorType >::is_non_negative | ( | ) | const |
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).
BlockVectorBase & BlockVectorBase< VectorType >::operator+= | ( | const BlockVectorBase< VectorType > & | V | ) |
Addition operator. Fast equivalent to U.add(1, V)
.
BlockVectorBase & BlockVectorBase< VectorType >::operator-= | ( | const BlockVectorBase< VectorType > & | V | ) |
Subtraction operator. Fast equivalent to U.add(-1, V)
.
void BlockVectorBase< VectorType >::add | ( | const std::vector< size_type > & | indices, |
const std::vector< Number > & | values | ||
) |
A collective add operation: This function adds a whole set of values stored in values
to the vector components specified by indices
.
void BlockVectorBase< VectorType >::add | ( | const std::vector< size_type > & | indices, |
const Vector< Number > & | values | ||
) |
This is a second collective add operation. As a difference, this function takes a deal.II vector of values.
void BlockVectorBase< VectorType >::add | ( | const size_type | n_elements, |
const size_type * | indices, | ||
const Number * | values | ||
) |
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.
void BlockVectorBase< VectorType >::add | ( | const value_type | s | ) |
\(U(0-DIM)+=s\). Addition of s
to all components. Note that s
is a scalar and not a vector.
void BlockVectorBase< VectorType >::add | ( | const value_type | a, |
const BlockVectorBase< VectorType > & | V | ||
) |
U+=a*V. Simple addition of a scaled vector.
void BlockVectorBase< VectorType >::add | ( | const value_type | a, |
const BlockVectorBase< VectorType > & | V, | ||
const value_type | b, | ||
const BlockVectorBase< VectorType > & | W | ||
) |
U+=a*V+b*W. Multiple addition of scaled vectors.
void BlockVectorBase< VectorType >::sadd | ( | const value_type | s, |
const BlockVectorBase< VectorType > & | V | ||
) |
U=s*U+V. Scaling and simple vector addition.
void BlockVectorBase< VectorType >::sadd | ( | const value_type | s, |
const value_type | a, | ||
const BlockVectorBase< VectorType > & | V | ||
) |
U=s*U+a*V. Scaling and simple addition.
void BlockVectorBase< VectorType >::sadd | ( | const value_type | s, |
const value_type | a, | ||
const BlockVectorBase< VectorType > & | V, | ||
const value_type | b, | ||
const BlockVectorBase< VectorType > & | W | ||
) |
U=s*U+a*V+b*W. Scaling and multiple addition.
void BlockVectorBase< VectorType >::sadd | ( | const value_type | s, |
const value_type | a, | ||
const BlockVectorBase< VectorType > & | V, | ||
const value_type | b, | ||
const BlockVectorBase< VectorType > & | W, | ||
const value_type | c, | ||
const BlockVectorBase< VectorType > & | X | ||
) |
U=s*U+a*V+b*W+c*X. Scaling and multiple addition.
BlockVectorBase & BlockVectorBase< VectorType >::operator*= | ( | const value_type | factor | ) |
Scale each element of the vector by a constant value.
BlockVectorBase & BlockVectorBase< VectorType >::operator/= | ( | const value_type | factor | ) |
Scale each element of the vector by the inverse of the given value.
void BlockVectorBase< VectorType >::scale | ( | const BlockVector2 & | v | ) |
Multiply each element of this vector by the corresponding element of v
.
void BlockVectorBase< VectorType >::equ | ( | const value_type | a, |
const BlockVector2 & | V | ||
) |
U=a*V. Assignment.
void BlockVectorBase< VectorType >::update_ghost_values | ( | ) | const |
Update the ghost values by calling update_ghost_values
for each block.
MPI_Comm BlockVectorBase< VectorType >::get_mpi_communicator | ( | ) | const |
This function returns the MPI communicator of the vector in the underlying blocks or, if the vector has not been initialized, the empty MPI_COMM_SELF.
std::size_t BlockVectorBase< VectorType >::memory_consumption | ( | ) | const |
Determine an estimate for the memory consumption (in bytes) of this object.
|
pure virtualinherited |
Extract a subset of the vector specified by indices
into the output array elements
.
|
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 131 of file enable_observer_pointer.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 151 of file enable_observer_pointer.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 322 of file enable_observer_pointer.h.
|
inlineinherited |
List the subscribers to the input stream
.
Definition at line 339 of file enable_observer_pointer.h.
|
inherited |
List the subscribers to deallog
.
Definition at line 199 of file enable_observer_pointer.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 331 of file enable_observer_pointer.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 53 of file enable_observer_pointer.cc.
|
friend |
Definition at line 979 of file block_vector_base.h.
Definition at line 982 of file block_vector_base.h.
|
protected |
Pointer to the array of components.
Definition at line 969 of file block_vector_base.h.
|
protected |
Object managing the transformation between global indices and indices within the different blocks.
Definition at line 975 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 227 of file enable_observer_pointer.h.
|
mutableprivateinherited |
In this map, we count subscriptions for each different identification string supplied to subscribe().
Definition at line 233 of file enable_observer_pointer.h.
|
mutableprivateinherited |
In this vector, we store pointers to the validity bool in the ObserverPointer objects that subscribe to this class.
Definition at line 249 of file enable_observer_pointer.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 257 of file enable_observer_pointer.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 280 of file enable_observer_pointer.h.