Reference documentation for deal.II version 9.3.3
|
#include <deal.II/lac/trilinos_vector.h>
Public Types | |
using | value_type = TrilinosScalar |
using | real_type = TrilinosScalar |
using | size_type = ::types::global_dof_index |
using | iterator = value_type * |
using | const_iterator = const value_type * |
using | reference = internal::VectorReference |
using | const_reference = const internal::VectorReference |
Public Member Functions | |
1: Basic Object-handling | |
Vector () | |
Vector (const Vector &v) | |
Vector (const IndexSet ¶llel_partitioning, const MPI_Comm &communicator=MPI_COMM_WORLD) | |
Vector (const IndexSet &local, const IndexSet &ghost, const MPI_Comm &communicator=MPI_COMM_WORLD) | |
Vector (const IndexSet ¶llel_partitioning, const Vector &v, const MPI_Comm &communicator=MPI_COMM_WORLD) | |
template<typename Number > | |
Vector (const IndexSet ¶llel_partitioning, const ::Vector< Number > &v, const MPI_Comm &communicator=MPI_COMM_WORLD) | |
Vector (Vector &&v) noexcept | |
~Vector () override=default | |
void | clear () |
void | reinit (const Vector &v, const bool omit_zeroing_entries=false, const bool allow_different_maps=false) |
void | reinit (const IndexSet ¶llel_partitioning, const MPI_Comm &communicator=MPI_COMM_WORLD, const bool omit_zeroing_entries=false) |
void | reinit (const IndexSet &locally_owned_entries, const IndexSet &ghost_entries, const MPI_Comm &communicator=MPI_COMM_WORLD, const bool vector_writable=false) |
void | reinit (const BlockVector &v, const bool import_data=false) |
void | compress (::VectorOperation::values operation) |
Vector & | operator= (const TrilinosScalar s) |
Vector & | operator= (const Vector &v) |
Vector & | operator= (Vector &&v) noexcept |
template<typename Number > | |
Vector & | operator= (const ::Vector< Number > &v) |
void | import_nonlocal_data_for_fe (const ::TrilinosWrappers::SparseMatrix &matrix, const Vector &vector) |
void | import (const LinearAlgebra::ReadWriteVector< double > &rwv, const VectorOperation::values operation) |
bool | operator== (const Vector &v) const |
bool | operator!= (const Vector &v) const |
size_type | size () const |
size_type | local_size () const |
size_type | locally_owned_size () const |
std::pair< size_type, size_type > | local_range () const |
bool | in_local_range (const size_type index) const |
IndexSet | locally_owned_elements () const |
bool | has_ghost_elements () const |
void | update_ghost_values () const |
TrilinosScalar | operator* (const Vector &vec) const |
real_type | norm_sqr () const |
TrilinosScalar | mean_value () const |
TrilinosScalar | min () const |
TrilinosScalar | max () const |
real_type | l1_norm () const |
real_type | l2_norm () const |
real_type | lp_norm (const TrilinosScalar p) const |
real_type | linfty_norm () const |
TrilinosScalar | add_and_dot (const TrilinosScalar a, const Vector &V, const Vector &W) |
bool | all_zero () const |
bool | is_non_negative () const |
2: Data-Access | |
reference | operator() (const size_type index) |
TrilinosScalar | operator() (const size_type index) const |
reference | operator[] (const size_type index) |
TrilinosScalar | operator[] (const size_type index) const |
void | extract_subvector_to (const std::vector< size_type > &indices, std::vector< TrilinosScalar > &values) const |
template<typename ForwardIterator , typename OutputIterator > | |
void | extract_subvector_to (ForwardIterator indices_begin, const ForwardIterator indices_end, OutputIterator values_begin) const |
iterator | begin () |
const_iterator | begin () const |
iterator | end () |
const_iterator | end () const |
3: Modification of vectors | |
void | set (const std::vector< size_type > &indices, const std::vector< TrilinosScalar > &values) |
void | set (const std::vector< size_type > &indices, const ::Vector< TrilinosScalar > &values) |
void | set (const size_type n_elements, const size_type *indices, const TrilinosScalar *values) |
void | add (const std::vector< size_type > &indices, const std::vector< TrilinosScalar > &values) |
void | add (const std::vector< size_type > &indices, const ::Vector< TrilinosScalar > &values) |
void | add (const size_type n_elements, const size_type *indices, const TrilinosScalar *values) |
Vector & | operator*= (const TrilinosScalar factor) |
Vector & | operator/= (const TrilinosScalar factor) |
Vector & | operator+= (const Vector &V) |
Vector & | operator-= (const Vector &V) |
void | add (const TrilinosScalar s) |
void | add (const Vector &V, const bool allow_different_maps=false) |
void | add (const TrilinosScalar a, const Vector &V) |
void | add (const TrilinosScalar a, const Vector &V, const TrilinosScalar b, const Vector &W) |
void | sadd (const TrilinosScalar s, const Vector &V) |
void | sadd (const TrilinosScalar s, const TrilinosScalar a, const Vector &V) |
void | scale (const Vector &scaling_factors) |
void | equ (const TrilinosScalar a, const Vector &V) |
Related Functions | |
(Note that these are not member functions.) | |
void | swap (Vector &u, Vector &v) |
4: Mixed stuff | |
Epetra_CombineMode | last_action |
bool | compressed |
bool | has_ghosts |
std::unique_ptr< Epetra_FEVector > | vector |
std::unique_ptr< Epetra_MultiVector > | nonlocal_vector |
IndexSet | owned_elements |
class | internal::VectorReference |
const Epetra_MultiVector & | trilinos_vector () const |
Epetra_FEVector & | trilinos_vector () |
const Epetra_BlockMap & | trilinos_partitioner () const |
void | print (std::ostream &out, const unsigned int precision=3, const bool scientific=true, const bool across=true) const |
void | swap (Vector &v) |
std::size_t | memory_consumption () const |
const MPI_Comm & | get_mpi_communicator () const |
static ::ExceptionBase & | ExcDifferentParallelPartitioning () |
static ::ExceptionBase & | ExcTrilinosError (int arg1) |
static ::ExceptionBase & | ExcAccessToNonLocalElement (size_type arg1, size_type arg2, size_type arg3, size_type arg4) |
Subscriptor functionality | |
Classes derived from Subscriptor provide a facility to subscribe to this object. This is mostly used by the SmartPointer class. | |
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 |
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 |
static std::mutex | mutex |
void | check_no_subscribers () const noexcept |
This class implements a wrapper to use the Trilinos distributed vector class Epetra_FEVector, the (parallel) partitioning of which is governed by an Epetra_Map. The Epetra_FEVector is precisely the kind of vector we deal with all the time - we probably get it from some assembly process, where also entries not locally owned might need to written and hence need to be forwarded to the owner.
The interface of this class is modeled after the existing Vector class in deal.II. It has almost the same member functions, and is often exchangeable. However, since Trilinos only supports a single scalar type (double), it is not templated, and only works with that type.
Note that Trilinos only guarantees that operations do what you expect if the function GlobalAssemble
has been called after vector assembly in order to distribute the data. This is necessary since some processes might have accumulated data of elements that are not owned by themselves, but must be sent to the owning process. In order to avoid using the wrong data, you need to call Vector::compress() before you actually use the vectors.
The parallel functionality of Trilinos is built on top of the Message Passing Interface (MPI). MPI's communication model is built on collective communications: if one process wants something from another, that other process has to be willing to accept this communication. A process cannot query data from another process by calling a remote function, without that other process expecting such a transaction. The consequence is that most of the operations in the base class of this class have to be called collectively. For example, if you want to compute the l2 norm of a parallel vector, all processes across which this vector is shared have to call the l2_norm
function. If you don't do this, but instead only call the l2_norm
function on one process, then the following happens: This one process will call one of the collective MPI functions and wait for all the other processes to join in on this. Since the other processes don't call this function, you will either get a time-out on the first process, or, worse, by the time the next a call to a Trilinos function generates an MPI message on the other processes, you will get a cryptic message that only a subset of processes attempted a communication. These bugs can be very hard to figure out, unless you are well-acquainted with the communication model of MPI, and know which functions may generate MPI messages.
One particular case, where an MPI message may be generated unexpectedly is discussed below.
Trilinos does of course allow read access to individual elements of a vector, but in the distributed case only to elements that are stored locally. We implement this through calls like d=vec(i)
. However, if you access an element outside the locally stored range, an exception is generated.
In contrast to read access, Trilinos (and the respective deal.II wrapper classes) allow to write (or add) to individual elements of vectors, even if they are stored on a different process. You can do this by writing into or adding to elements using the syntax vec(i)=d
or vec(i)+=d
, or similar operations. There is one catch, however, that may lead to very confusing error messages: Trilinos requires application programs to call the compress() function when they switch from performing a set of operations that add to elements, to performing a set of operations that write to elements. The reasoning is that all processes might accumulate addition operations to elements, even if multiple processes write to the same elements. By the time we call compress() the next time, all these additions are executed. However, if one process adds to an element, and another overwrites to it, the order of execution would yield non-deterministic behavior if we don't make sure that a synchronization with compress() happens in between.
In order to make sure these calls to compress() happen at the appropriate time, the deal.II wrappers keep a state variable that store which is the presently allowed operation: additions or writes. If it encounters an operation of the opposite kind, it calls compress() and flips the state. This can sometimes lead to very confusing behavior, in code that may for example look like this:
This code can run into trouble: by the time we see the first addition operation, we need to flush the overwrite buffers for the vector, and the deal.II library will do so by calling compress(). However, it will only do so for all processes that actually do an addition – if the condition is never true for one of the processes, then this one will not get to the actual compress() call, whereas all the other ones do. This gets us into trouble, since all the other processes hang in the call to flush the write buffers, while the one other process advances to the call to compute the l2 norm. At this time, you will get an error that some operation was attempted by only a subset of processes. This behavior may seem surprising, unless you know that write/addition operations on single elements may trigger this behavior.
The problem described here may be avoided by placing additional calls to compress(), or making sure that all processes do the same type of operations at the same time, for example by placing zero additions if necessary.
Parallel vectors come in two kinds: without and with ghost elements. Vectors without ghost elements uniquely partition the vector elements between processors: each vector entry has exactly one processor that owns it. For such vectors, you can read those elements that the processor you are currently on owns, and you can write into any element whether you own it or not: if you don't own it, the value written or added to a vector element will be shipped to the processor that owns this vector element the next time you call compress(), as described above.
What we call a 'ghosted' vector (see vectors with ghost elements ) is simply a view of the parallel vector where the element distributions overlap. The 'ghosted' Trilinos vector in itself has no idea of which entries are ghosted and which are locally owned. In fact, a ghosted vector may not even store all of the elements a non- ghosted vector would store on the current processor. Consequently, for Trilinos vectors, there is no notion of an 'owner' of vector elements in the way we have it in the non-ghost case view.
This explains why we do not allow writing into ghosted vectors on the Trilinos side: Who would be responsible for taking care of the duplicated entries, given that there is not such information as locally owned indices? In other words, since a processor doesn't know which other processors own an element, who would it send a value to if one were to write to it? The only possibility would be to send this information to all other processors, but that is clearly not practical. Thus, we only allow reading from ghosted vectors, which however we do very often.
So how do you fill a ghosted vector if you can't write to it? This only happens through the assignment with a non-ghosted vector. It can go both ways (non-ghosted is assigned to a ghosted vector, and a ghosted vector is assigned to a non-ghosted one; the latter one typically only requires taking out the locally owned part as most often ghosted vectors store a superset of elements of non-ghosted ones). In general, you send data around with that operation and it all depends on the different views of the two vectors. Trilinos also allows you to get subvectors out of a big vector that way.
When writing into Trilinos vectors from several threads in shared memory, several things must be kept in mind as there is no built-in locks in this class to prevent data races. Simultaneous access to the same vector entry at the same time results in data races and must be explicitly avoided by the user. However, it is possible to access different entries of the vector from several threads simultaneously when only one MPI process is present or the vector has been constructed with an additional index set for ghost entries in write mode.
2008, 2009, 2017
Definition at line 403 of file trilinos_vector.h.