Reference documentation for deal.II version 9.2.0
\(\newcommand{\dealvcentcolon}{\mathrel{\mathop{:}}}\) \(\newcommand{\dealcoloneq}{\dealvcentcolon\mathrel{\mkern-1.2mu}=}\) \(\newcommand{\jump}[1]{\left[\!\left[ #1 \right]\!\right]}\) \(\newcommand{\average}[1]{\left\{\!\left\{ #1 \right\}\!\right\}}\)
Public Types | Public Member Functions | List of all members
LinearAlgebra::distributed::Vector< Number, MemorySpace > Class Template Reference

#include <deal.II/lac/la_parallel_vector.h>

Inheritance diagram for LinearAlgebra::distributed::Vector< Number, MemorySpace >:
[legend]

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 >
Vectoroperator= (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 ()
 
Subscriptoroperator= (const Subscriptor &)
 
Subscriptoroperator= (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)
 

6: Mixed stuff

std::shared_ptr< const Utilities::MPI::Partitionerpartitioner
 
size_type allocated_size
 
mutable ::MemorySpace::MemorySpaceData< Number, MemorySpace > data
 
std::shared_ptr<::parallel::internal::TBBPartitionerthread_loop_partitioner
 
mutable ::MemorySpace::MemorySpaceData< Number, MemorySpace > import_data
 
bool vector_is_ghosted
 
std::vector< MPI_Requestcompress_requests
 
std::vector< MPI_Requestupdate_ghost_values_requests
 
std::mutex mutex
 
template<typename Number2 , typename MemorySpace2 >
class Vector
 
template<typename Number2 >
class BlockVector
 
const MPI_Commget_mpi_communicator () const
 
const std::shared_ptr< const Utilities::MPI::Partitioner > & get_partitioner () const
 
bool partitioners_are_compatible (const Utilities::MPI::Partitioner &part) const
 
bool partitioners_are_globally_compatible (const Utilities::MPI::Partitioner &part) const
 
void set_ghost_state (const bool ghosted) const
 
static ::ExceptionBaseExcVectorTypeNotCompatible ()
 
static ::ExceptionBaseExcNotAllowedForCuda ()
 
static ::ExceptionBaseExcNonMatchingElements (Number arg1, Number arg2, unsigned int arg3)
 
static ::ExceptionBaseExcAccessToNonLocalElement (size_type arg1, size_type arg2, size_type arg3, size_type arg4)
 
void add_local (const Number a, const VectorSpaceVector< Number > &V)
 
void sadd_local (const Number s, const Number a, const VectorSpaceVector< Number > &V)
 
template<typename Number2 >
Number inner_product_local (const Vector< Number2, MemorySpace > &V) const
 
real_type norm_sqr_local () const
 
Number mean_value_local () const
 
real_type l1_norm_local () const
 
real_type lp_norm_local (const real_type p) const
 
real_type linfty_norm_local () const
 
Number add_and_dot_local (const Number a, const Vector< Number, MemorySpace > &V, const Vector< Number, MemorySpace > &W)
 
void clear_mpi_requests ()
 
void resize_val (const size_type new_allocated_size)
 

Additional Inherited Members

- Static Public Member Functions inherited from Subscriptor
static ::ExceptionBaseExcInUse (int arg1, std::string arg2, std::string arg3)
 
static ::ExceptionBaseExcNoSubscriber (std::string arg1, std::string arg2)
 

Detailed Description

template<typename Number, typename MemorySpace = MemorySpace::Host>
class LinearAlgebra::distributed::Vector< Number, MemorySpace >

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:

Functions related to parallel functionality:

This vector can take two different states with respect to ghost elements:

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.

Limitations regarding the vector size

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.

CUDA support

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.

Note
By default, all the ranks will try to access the device 0. This is fine is if you have one rank per node and one gpu per node. If you have multiple GPUs on one node, we need each process to access a different GPU. If each node has the same number of GPUs, this can be done as follows: int n_devices = 0; cudaGetDeviceCount(&n_devices); int device_id = my_rank % n_devices; cudaSetDevice(device_id);
See also
CUDAWrappers
Author
Katharina Kormann, Martin Kronbichler, Bruno Turcksin 2010, 2011, 2016, 2018

Definition at line 226 of file la_parallel_vector.h.

Member Function Documentation

◆ operator=()

template<typename Number , typename MemorySpace = MemorySpace::Host>
template<typename number >
Vector& LinearAlgebra::distributed::Vector< Number, MemorySpace >::operator= ( const ::Vector< number > &  v)

Definition at line 509 of file trilinos_vector.cc.


The documentation for this class was generated from the following files:
Utilities::CUDA::copy_to_dev
void copy_to_dev(const std::vector< T > &vector_host, T *pointer_dev)
Definition: cuda.h:146
Vector
Definition: mapping_q1_eulerian.h:32