Reference documentation for deal.II version 9.2.0
|
#include <deal.II/lac/read_write_vector.h>
Classes | |
class | FunctorTemplate |
Public Types | |
using | value_type = Number |
using | pointer = value_type * |
using | const_pointer = const value_type * |
using | iterator = value_type * |
using | const_iterator = const value_type * |
using | reference = value_type & |
using | const_reference = const value_type & |
using | size_type = types::global_dof_index |
using | real_type = typename numbers::NumberTraits< Number >::real_type |
Public Member Functions | |
1: Basic Object-handling | |
ReadWriteVector () | |
ReadWriteVector (const ReadWriteVector< Number > &in_vector) | |
ReadWriteVector (const size_type size) | |
ReadWriteVector (const IndexSet &locally_stored_indices) | |
~ReadWriteVector () override=default | |
virtual void | reinit (const size_type size, const bool omit_zeroing_entries=false) |
template<typename Number2 > | |
void | reinit (const ReadWriteVector< Number2 > &in_vector, const bool omit_zeroing_entries=false) |
virtual void | reinit (const IndexSet &locally_stored_indices, const bool omit_zeroing_entries=false) |
void | reinit (const TrilinosWrappers::MPI::Vector &trilinos_vec) |
template<typename Functor > | |
void | apply (const Functor &func) |
void | swap (ReadWriteVector< Number > &v) |
ReadWriteVector< Number > & | operator= (const ReadWriteVector< Number > &in_vector) |
template<typename Number2 > | |
ReadWriteVector< Number > & | operator= (const ReadWriteVector< Number2 > &in_vector) |
ReadWriteVector< Number > & | operator= (const Number s) |
template<typename MemorySpace > | |
void | import (const distributed::Vector< Number, MemorySpace > &vec, VectorOperation::values operation, const std::shared_ptr< const CommunicationPatternBase > &communication_pattern=std::shared_ptr< const CommunicationPatternBase >()) |
void | import (const PETScWrappers::MPI::Vector &petsc_vec, VectorOperation::values operation, const std::shared_ptr< const CommunicationPatternBase > &communication_pattern=std::shared_ptr< const CommunicationPatternBase >()) |
void | import (const TrilinosWrappers::MPI::Vector &trilinos_vec, VectorOperation::values operation, const std::shared_ptr< const CommunicationPatternBase > &communication_pattern=std::shared_ptr< const CommunicationPatternBase >()) |
void | import (const TpetraWrappers::Vector< Number > &tpetra_vec, VectorOperation::values operation, const std::shared_ptr< const CommunicationPatternBase > &communication_pattern=std::shared_ptr< const CommunicationPatternBase >()) |
void | import (const EpetraWrappers::Vector &epetra_vec, VectorOperation::values operation, const std::shared_ptr< const CommunicationPatternBase > &communication_pattern=std::shared_ptr< const CommunicationPatternBase >()) |
void | import (const CUDAWrappers::Vector< Number > &cuda_vec, VectorOperation::values operation, const std::shared_ptr< const CommunicationPatternBase > &communication_pattern=std::shared_ptr< const CommunicationPatternBase >()) |
size_type | size () const |
size_type | n_elements () const |
const IndexSet & | get_stored_elements () const |
iterator | begin () |
const_iterator | begin () const |
iterator | end () |
const_iterator | end () const |
2: Data-Access | |
Number | operator() (const size_type global_index) const |
Number & | operator() (const size_type global_index) |
Number | operator[] (const size_type global_index) const |
Number & | operator[] (const size_type global_index) |
template<typename Number2 > | |
void | extract_subvector_to (const std::vector< size_type > &indices, std::vector< Number2 > &values) const |
template<typename ForwardIterator , typename OutputIterator > | |
void | extract_subvector_to (ForwardIterator indices_begin, const ForwardIterator indices_end, OutputIterator values_begin) const |
Number | local_element (const size_type local_index) const |
Number & | local_element (const size_type local_index) |
Public Member Functions inherited from Subscriptor | |
Subscriptor () | |
Subscriptor (const Subscriptor &) | |
Subscriptor (Subscriptor &&) noexcept | |
virtual | ~Subscriptor () |
Subscriptor & | operator= (const Subscriptor &) |
Subscriptor & | operator= (Subscriptor &&) noexcept |
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) |
3: Modification of vectors | |
IndexSet | stored_elements |
IndexSet | source_stored_elements |
std::shared_ptr< CommunicationPatternBase > | comm_pattern |
std::unique_ptr< Number[], decltype(std::free) * > | values |
std::shared_ptr<::parallel::internal::TBBPartitioner > | thread_loop_partitioner |
template<typename Number2 > | |
class | ReadWriteVector |
template<typename Number2 > | |
void | add (const std::vector< size_type > &indices, const std::vector< Number2 > &values) |
template<typename Number2 > | |
void | add (const std::vector< size_type > &indices, const ReadWriteVector< Number2 > &values) |
template<typename Number2 > | |
void | add (const size_type n_elements, const size_type *indices, const Number2 *values) |
void | print (std::ostream &out, const unsigned int precision=3, const bool scientific=true) const |
std::size_t | memory_consumption () const |
void | import (const Tpetra::Vector< Number, int, types::global_dof_index > &tpetra_vector, const IndexSet &locally_owned_elements, VectorOperation::values operation, const MPI_Comm &mpi_comm, const std::shared_ptr< const CommunicationPatternBase > &communication_pattern) |
void | import (const Epetra_MultiVector &multivector, const IndexSet &locally_owned_elements, VectorOperation::values operation, const MPI_Comm &mpi_comm, const std::shared_ptr< const CommunicationPatternBase > &communication_pattern) |
unsigned int | global_to_local (const types::global_dof_index global_index) const |
void | resize_val (const size_type new_allocated_size) |
TpetraWrappers::CommunicationPattern | create_tpetra_comm_pattern (const IndexSet &source_index_set, const MPI_Comm &mpi_comm) |
EpetraWrappers::CommunicationPattern | create_epetra_comm_pattern (const IndexSet &source_index_set, const MPI_Comm &mpi_comm) |
Additional Inherited Members | |
Static Public Member Functions inherited from Subscriptor | |
static ::ExceptionBase & | ExcInUse (int arg1, std::string arg2, std::string arg3) |
static ::ExceptionBase & | ExcNoSubscriber (std::string arg1, std::string arg2) |
ReadWriteVector is intended to represent vectors in \({\mathbb R}^N\) for which it stores all or a subset of elements. The latter case in important in parallel computations, where \(N\) may be so large that no processor can actually all elements of a solution vector, but where this is also not necessary: one typically only has to store the values of degrees of freedom that live on cells that are locally owned plus potentially those degrees of freedom that live on ghost cells.
This class allows to access individual elements to be read or written. However, it does not allow global operations such as taking the norm. ReadWriteVector can be used to read and write elements in vectors derived from VectorSpaceVector such as TrilinosWrappers::MPI::Vector and PETScWrappers::MPI::Vector.
Most of the time, one will simply read from or write into a vector of the current class using the global numbers of these degrees of freedom. This is done using operator()() or operator[]() which call global_to_local() to transform the global index into a local one. In such cases, it is clear that one can only access elements of the vector that the current object indeed stores.
However, it is also possible to access elements in the order in which they are stored by the current object. In other words, one is not interested in accessing elements with their global indices, but instead using an enumeration that only takes into account the elements that are actually stored. This is facilitated by the local_element() function. To this end, it is necessary to know in which order the current class stores its element. The elements of all the consecutive ranges are stored in ascending order of the first index of each range. The function IndexSet::largest_range_starting_index() can be used to get the first index of the largest range.
Definition at line 131 of file read_write_vector.h.