Reference documentation for deal.II version 9.2.0
|
#include <deal.II/lac/la_parallel_vector.h>
Public Types | |
using | memory_space = MemorySpace |
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 Types inherited from LinearAlgebra::VectorSpaceVector< Number > | |
using | value_type = Number |
using | size_type = types::global_dof_index |
using | real_type = typename numbers::NumberTraits< Number >::real_type |
Public Member Functions | |
template<typename number > | |
Vector & | operator= (const ::Vector< number > &v) |
1: Basic Object-handling | |
Vector () | |
Vector (const Vector< Number, MemorySpace > &in_vector) | |
Vector (const size_type size) | |
Vector (const IndexSet &local_range, const IndexSet &ghost_indices, const MPI_Comm communicator) | |
Vector (const IndexSet &local_range, const MPI_Comm communicator) | |
Vector (const std::shared_ptr< const Utilities::MPI::Partitioner > &partitioner) | |
virtual | ~Vector () override |
void | reinit (const size_type size, const bool omit_zeroing_entries=false) |
template<typename Number2 > | |
void | reinit (const Vector< Number2, MemorySpace > &in_vector, const bool omit_zeroing_entries=false) |
void | reinit (const IndexSet &local_range, const IndexSet &ghost_indices, const MPI_Comm communicator) |
void | reinit (const IndexSet &local_range, const MPI_Comm communicator) |
void | reinit (const std::shared_ptr< const Utilities::MPI::Partitioner > &partitioner) |
void | swap (Vector< Number, MemorySpace > &v) |
Vector< Number, MemorySpace > & | operator= (const Vector< Number, MemorySpace > &in_vector) |
template<typename Number2 > | |
Vector< Number, MemorySpace > & | operator= (const Vector< Number2, MemorySpace > &in_vector) |
2: Parallel data exchange | |
virtual void | compress (::VectorOperation::values operation) override |
void | update_ghost_values () const |
void | compress_start (const unsigned int communication_channel=0, ::VectorOperation::values operation=VectorOperation::add) |
void | compress_finish (::VectorOperation::values operation) |
void | update_ghost_values_start (const unsigned int communication_channel=0) const |
void | update_ghost_values_finish () const |
void | zero_out_ghosts () const |
bool | has_ghost_elements () const |
template<typename Number2 > | |
void | copy_locally_owned_data_from (const Vector< Number2, MemorySpace > &src) |
template<typename MemorySpace2 > | |
void | import (const Vector< Number, MemorySpace2 > &src, VectorOperation::values operation) |
3: Implementation of VectorSpaceVector | |
virtual void | reinit (const VectorSpaceVector< Number > &V, const bool omit_zeroing_entries=false) override |
virtual Vector< Number, MemorySpace > & | operator*= (const Number factor) override |
virtual Vector< Number, MemorySpace > & | operator/= (const Number factor) override |
virtual Vector< Number, MemorySpace > & | operator+= (const VectorSpaceVector< Number > &V) override |
virtual Vector< Number, MemorySpace > & | operator-= (const VectorSpaceVector< Number > &V) override |
virtual void | import (const LinearAlgebra::ReadWriteVector< Number > &V, VectorOperation::values operation, std::shared_ptr< const CommunicationPatternBase > communication_pattern=std::shared_ptr< const CommunicationPatternBase >()) override |
virtual Number | operator* (const VectorSpaceVector< Number > &V) const override |
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 |
4: Other vector operations not included in VectorSpaceVector | |
virtual Vector< Number, MemorySpace > & | operator= (const Number s) override |
template<typename OtherNumber > | |
void | add (const std::vector< size_type > &indices, const ::Vector< OtherNumber > &values) |
template<typename OtherNumber > | |
void | add (const size_type n_elements, const size_type *indices, const OtherNumber *values) |
void | sadd (const Number s, const Vector< Number, MemorySpace > &V) |
5: Entry access and local data representation | |
size_type | local_size () const |
bool | in_local_range (const size_type global_index) const |
iterator | begin () |
const_iterator | begin () const |
iterator | end () |
const_iterator | end () const |
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) |
Number | local_element (const size_type local_index) const |
Number & | local_element (const size_type local_index) |
Number * | get_values () const |
template<typename OtherNumber > | |
void | extract_subvector_to (const std::vector< size_type > &indices, std::vector< OtherNumber > &values) const |
template<typename ForwardIterator , typename OutputIterator > | |
void | extract_subvector_to (ForwardIterator indices_begin, const ForwardIterator indices_end, OutputIterator values_begin) const |
virtual bool | all_zero () const override |
virtual Number | mean_value () const override |
real_type | lp_norm (const real_type p) const |
Public Member Functions inherited from LinearAlgebra::VectorSpaceVector< Number > | |
virtual void | compress (VectorOperation::values) |
virtual | ~VectorSpaceVector ()=default |
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) |
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) |
Implementation of a parallel vector class. The design of this class is similar to the standard Vector class in deal.II, with the exception that storage is distributed with MPI.
The vector is designed for the following scheme of parallel partitioning:
[my_first_index,my_last_index)
. reinit (locally_owned, ghost_indices, communicator)
, and retained until the partitioning is changed. This allows for efficient parallel communication of indices. In particular, it stores the communication pattern, rather than having to compute it again for every communication. For more information on ghost vectors, see also the glossary entry on vectors with ghost elements. local_element()
. Locally owned indices are placed first, [0, local_size()), and then all ghost indices follow after them contiguously, [local_size(), local_size()+n_ghost_entries()). Functions related to parallel functionality:
compress()
goes through the data associated with ghost indices and communicates it to the owner process, which can then add it to the correct position. This can be used e.g. after having run an assembly routine involving ghosts that fill this vector. Note that the insert
mode of compress()
does not set the elements included in ghost entries but simply discards them, assuming that the owning processor has set them to the desired value already (See also the glossary entry on compress). update_ghost_values()
function imports the data from the owning processor to the ghost indices in order to provide read access to the data associated with ghosts. This vector can take two different states with respect to ghost elements:
operator= (0.)
), the vector does only allow writing into ghost elements but not reading from ghost elements. true
exactly when ghost elements have been updated and false
otherwise, irrespective of the actual number of ghost entries in the vector layout (for that information, use n_ghost_entries() instead). This vector uses the facilities of the class Vector<Number> for implementing the operations on the local range of the vector. In particular, it also inherits thread parallelism that splits most vector-vector operations into smaller chunks if the program uses multiple threads. This may or may not be desired when working also with MPI.
This vector class is based on two different number types for indexing. The so-called global index type encodes the overall size of the vector. Its type is types::global_dof_index. The largest possible value is 2^32-1
or approximately 4 billion in case 64 bit integers are disabled at configuration of deal.II (default case) or 2^64-1
or approximately 10^19
if 64 bit integers are enabled (see the glossary entry on When to use types::global_dof_index instead of unsigned int for further information).
The second relevant index type is the local index used within one MPI rank. As opposed to the global index, the implementation assumes 32-bit unsigned integers unconditionally. In other words, to actually use a vector with more than four billion entries, you need to use MPI with more than one rank (which in general is a safe assumption since four billion entries consume at least 16 GB of memory for floats or 32 GB of memory for doubles) and enable 64-bit indices. If more than 4 billion local elements are present, the implementation tries to detect that, which triggers an exception and aborts the code. Note, however, that the detection of overflow is tricky and the detection mechanism might fail in some circumstances. Therefore, it is strongly recommended to not rely on this class to automatically detect the unsupported case.
This vector class supports two different memory spaces: Host and CUDA. By default, the memory space is Host and all the data are allocated on the CPU. When the memory space is CUDA, all the data is allocated on the GPU. The operations on the vector are performed on the chosen memory space. * From the host, there are two methods to access the elements of the Vector when using the CUDA memory space:
The import method is a lot safer and will perform an MPI communication if necessary. Since an MPI communication may be performed, import needs to be called on all the processors.
int n_devices = 0; cudaGetDeviceCount(&n_devices); int device_id = my_rank % n_devices; cudaSetDevice(device_id);
Definition at line 226 of file la_parallel_vector.h.
Vector& LinearAlgebra::distributed::Vector< Number, MemorySpace >::operator= | ( | const ::Vector< number > & | v | ) |
Definition at line 509 of file trilinos_vector.cc.