Reference documentation for deal.II version GIT relicensing-245-g36f19064f7 2024-03-29 07:20:02+00:00
\(\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\}}\)
Loading...
Searching...
No Matches
Classes | Public Types | Public Member Functions | Static Public Member Functions | Protected Member Functions | Protected Attributes | Private Types | Private Member Functions | Private Attributes | Static Private Attributes | Friends | List of all members
LinearAlgebra::ReadWriteVector< Number > Class Template Reference

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

Inheritance diagram for LinearAlgebra::ReadWriteVector< Number >:
Inheritance graph
[legend]

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

template<class Archive >
void serialize (Archive &ar, const unsigned int version)
 
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)
 
void import_elements (const ::Vector< Number > &vec, VectorOperation::values operation, const std::shared_ptr< const Utilities::MPI::CommunicationPatternBase > &communication_pattern={})
 
void import (const ::Vector< Number > &V, VectorOperation::values operation, const std::shared_ptr< const Utilities::MPI::CommunicationPatternBase > &communication_pattern={})
 
template<typename MemorySpace >
void import_elements (const distributed::Vector< Number, MemorySpace > &vec, VectorOperation::values operation, const std::shared_ptr< const Utilities::MPI::CommunicationPatternBase > &communication_pattern={})
 
template<typename MemorySpace >
void import (const distributed::Vector< Number, MemorySpace > &V, VectorOperation::values operation, const std::shared_ptr< const Utilities::MPI::CommunicationPatternBase > &communication_pattern={})
 
void import_elements (const PETScWrappers::MPI::Vector &petsc_vec, VectorOperation::values operation, const std::shared_ptr< const Utilities::MPI::CommunicationPatternBase > &communication_pattern={})
 
void import (const PETScWrappers::MPI::Vector &V, VectorOperation::values operation, const std::shared_ptr< const Utilities::MPI::CommunicationPatternBase > &communication_pattern={})
 
void import_elements (const TrilinosWrappers::MPI::Vector &trilinos_vec, VectorOperation::values operation, const std::shared_ptr< const Utilities::MPI::CommunicationPatternBase > &communication_pattern={})
 
void import (const TrilinosWrappers::MPI::Vector &V, VectorOperation::values operation, const std::shared_ptr< const Utilities::MPI::CommunicationPatternBase > &communication_pattern={})
 
template<typename Dummy = Number>
std::enable_if_t< std::is_same_v< Dummy, Number > &&::is_tpetra_type< Number >::value > import_elements (const TpetraWrappers::Vector< Number > &tpetra_vec, VectorOperation::values operation, const std::shared_ptr< const Utilities::MPI::CommunicationPatternBase > &communication_pattern={})
 
template<typename Dummy = Number>
std::enable_if_t< std::is_same_v< Dummy, Number > &&::is_tpetra_type< Number >::value > import (const TpetraWrappers::Vector< Number > &V, VectorOperation::values operation, const std::shared_ptr< const Utilities::MPI::CommunicationPatternBase > &communication_pattern={})
 
void import_elements (const EpetraWrappers::Vector &epetra_vec, VectorOperation::values operation, const std::shared_ptr< const Utilities::MPI::CommunicationPatternBase > &communication_pattern={})
 
void import (const EpetraWrappers::Vector &V, VectorOperation::values operation, const std::shared_ptr< const Utilities::MPI::CommunicationPatternBase > &communication_pattern={})
 
void import_elements (const CUDAWrappers::Vector< Number > &cuda_vec, VectorOperation::values operation, const std::shared_ptr< const Utilities::MPI::CommunicationPatternBase > &communication_pattern={})
 
void import (const CUDAWrappers::Vector< Number > &V, VectorOperation::values operation, const std::shared_ptr< const Utilities::MPI::CommunicationPatternBase > &communication_pattern={})
 
size_type size () const override
 
size_type locally_owned_size () const
 
const IndexSetget_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
 
virtual void extract_subvector_to (const ArrayView< const types::global_dof_index > &indices, ArrayView< Number > &entries) const override
 
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)
 
3: Modification of vectors
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
 
Subscriptor functionality

Classes derived from Subscriptor provide a facility to subscribe to this object. This is mostly used by the SmartPointer 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 ::ExceptionBaseExcInUse (int arg1, std::string arg2, std::string arg3)
 
static ::ExceptionBaseExcNoSubscriber (std::string arg1, std::string arg2)
 

Protected Member Functions

template<typename Dummy = Number>
std::enable_if_t< std::is_same_v< Dummy, Number > &&::is_tpetra_type< Number >::value > import_elements (const Tpetra::Vector< Number, int, types::signed_global_dof_index > &tpetra_vector, const IndexSet &locally_owned_elements, VectorOperation::values operation, const MPI_Comm mpi_comm, const std::shared_ptr< const Utilities::MPI::CommunicationPatternBase > &communication_pattern)
 
void import_elements (const Epetra_MultiVector &multivector, const IndexSet &locally_owned_elements, VectorOperation::values operation, const MPI_Comm mpi_comm, const std::shared_ptr< const Utilities::MPI::CommunicationPatternBase > &communication_pattern)
 
unsigned int global_to_local (const types::global_dof_index global_index) const
 
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)
 

Protected Attributes

IndexSet stored_elements
 
IndexSet source_stored_elements
 
std::shared_ptr< Utilities::MPI::CommunicationPatternBasecomm_pattern
 
AlignedVector< Number > values
 
std::shared_ptr<::parallel::internal::TBBPartitionerthread_loop_partitioner
 

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 intcounter
 
std::map< std::string, unsigned intcounter_map
 
std::vector< std::atomic< bool > * > validity_pointers
 
const std::type_info * object_info
 

Static Private Attributes

static std::mutex mutex
 

Friends

template<typename Number2 >
class ReadWriteVector
 

Detailed Description

template<typename Number>
class LinearAlgebra::ReadWriteVector< Number >

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 store 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 provides access to individual elements to be read or written. However, it does not allow global operations such as taking the norm or dot products between vectors.

Storing elements

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 136 of file read_write_vector.h.

Member Typedef Documentation

◆ value_type

template<typename Number >
using LinearAlgebra::ReadWriteVector< Number >::value_type = Number

Declare standard types used in all containers. These types parallel those in the C++ standard libraries vector<...> class.

Definition at line 144 of file read_write_vector.h.

◆ pointer

template<typename Number >
using LinearAlgebra::ReadWriteVector< Number >::pointer = value_type *

Definition at line 145 of file read_write_vector.h.

◆ const_pointer

template<typename Number >
using LinearAlgebra::ReadWriteVector< Number >::const_pointer = const value_type *

Definition at line 146 of file read_write_vector.h.

◆ iterator

template<typename Number >
using LinearAlgebra::ReadWriteVector< Number >::iterator = value_type *

Definition at line 147 of file read_write_vector.h.

◆ const_iterator

template<typename Number >
using LinearAlgebra::ReadWriteVector< Number >::const_iterator = const value_type *

Definition at line 148 of file read_write_vector.h.

◆ reference

template<typename Number >
using LinearAlgebra::ReadWriteVector< Number >::reference = value_type &

Definition at line 149 of file read_write_vector.h.

◆ const_reference

template<typename Number >
using LinearAlgebra::ReadWriteVector< Number >::const_reference = const value_type &

Definition at line 150 of file read_write_vector.h.

◆ size_type

template<typename Number >
using LinearAlgebra::ReadWriteVector< Number >::size_type = types::global_dof_index

Definition at line 151 of file read_write_vector.h.

◆ real_type

template<typename Number >
using LinearAlgebra::ReadWriteVector< Number >::real_type = typename numbers::NumberTraits<Number>::real_type

Definition at line 152 of file read_write_vector.h.

◆ map_value_type

using Subscriptor::map_value_type = decltype(counter_map)::value_type
privateinherited

The data type used in counter_map.

Definition at line 229 of file subscriptor.h.

◆ map_iterator

using Subscriptor::map_iterator = decltype(counter_map)::iterator
privateinherited

The iterator type used in counter_map.

Definition at line 234 of file subscriptor.h.

Constructor & Destructor Documentation

◆ ReadWriteVector() [1/4]

template<typename Number >
LinearAlgebra::ReadWriteVector< Number >::ReadWriteVector ( )

Empty constructor.

◆ ReadWriteVector() [2/4]

template<typename Number >
LinearAlgebra::ReadWriteVector< Number >::ReadWriteVector ( const ReadWriteVector< Number > &  in_vector)

Copy constructor.

◆ ReadWriteVector() [3/4]

template<typename Number >
LinearAlgebra::ReadWriteVector< Number >::ReadWriteVector ( const size_type  size)
explicit

Construct a vector given the size, the stored elements have their index in [0,size).

◆ ReadWriteVector() [4/4]

template<typename Number >
LinearAlgebra::ReadWriteVector< Number >::ReadWriteVector ( const IndexSet locally_stored_indices)
explicit

Construct a vector whose stored elements indices are given by the IndexSet locally_stored_indices.

◆ ~ReadWriteVector()

template<typename Number >
LinearAlgebra::ReadWriteVector< Number >::~ReadWriteVector ( )
overridedefault

Destructor.

Member Function Documentation

◆ reinit() [1/4]

template<typename Number >
virtual void LinearAlgebra::ReadWriteVector< Number >::reinit ( const size_type  size,
const bool  omit_zeroing_entries = false 
)
virtual

Set the global size of the vector to size. The stored elements have their index in [0,size).

If the flag omit_zeroing_entries is set to false, the memory will be initialized with zero, otherwise the memory will be untouched (and the user must make sure to fill it with reasonable data before using it).

◆ reinit() [2/4]

template<typename Number >
template<typename Number2 >
void LinearAlgebra::ReadWriteVector< Number >::reinit ( const ReadWriteVector< Number2 > &  in_vector,
const bool  omit_zeroing_entries = false 
)

Uses the same IndexSet as the one of the input vector in_vector and allocates memory for this vector.

If the flag omit_zeroing_entries is set to false, the memory will be initialized with zero, otherwise the memory will be untouched (and the user must make sure to fill it with reasonable data before using it).

◆ reinit() [3/4]

template<typename Number >
virtual void LinearAlgebra::ReadWriteVector< Number >::reinit ( const IndexSet locally_stored_indices,
const bool  omit_zeroing_entries = false 
)
virtual

Initializes the vector. The indices are specified by locally_stored_indices.

If the flag omit_zeroing_entries is set to false, the memory will be initialized with zero, otherwise the memory will be untouched (and the user must make sure to fill it with reasonable data before using it). locally_stored_indices.

◆ reinit() [4/4]

template<typename Number >
void LinearAlgebra::ReadWriteVector< Number >::reinit ( const TrilinosWrappers::MPI::Vector trilinos_vec)

Initialize this ReadWriteVector by supplying access to all locally available entries in the given ghosted or non-ghosted vector.

Note
This function currently copies the values from the argument into the ReadWriteVector, so modifications here will not modify trilinos_vec.

This function is mainly written for backwards-compatibility to get element access to a ghosted TrilinosWrappers::MPI::Vector inside the library.

◆ apply()

template<typename Number >
template<typename Functor >
void LinearAlgebra::ReadWriteVector< Number >::apply ( const Functor &  func)

Apply the functor func to each element of the vector. The functor should look like

struct Functor
{
void operator() (Number &value);
};
Note
This function requires that the header read_write_vector.templates.h be included.

◆ swap()

template<typename Number >
void LinearAlgebra::ReadWriteVector< Number >::swap ( ReadWriteVector< Number > &  v)

Swap the contents of this vector and the other vector v. One could do this operation with a temporary variable and copying over the data elements, but this function is significantly more efficient since it only swaps the pointers to the data of the two vectors and therefore does not need to allocate temporary storage and move data around.

This function is analogous to the swap function of all C++ standard containers. Also, there is a global function swap(u,v) that simply calls u.swap(v), again in analogy to standard functions.

◆ operator=() [1/3]

template<typename Number >
ReadWriteVector< Number > & LinearAlgebra::ReadWriteVector< Number >::operator= ( const ReadWriteVector< Number > &  in_vector)

Copies the data and the IndexSet of the input vector in_vector.

◆ operator=() [2/3]

template<typename Number >
template<typename Number2 >
ReadWriteVector< Number > & LinearAlgebra::ReadWriteVector< Number >::operator= ( const ReadWriteVector< Number2 > &  in_vector)

Copies the data and the IndexSet of the input vector in_vector.

◆ operator=() [3/3]

template<typename Number >
ReadWriteVector< Number > & LinearAlgebra::ReadWriteVector< Number >::operator= ( const Number  s)

Sets all elements of the vector to the scalar s. This operation is only allowed if s is equal to zero.

◆ import_elements() [1/9]

template<typename Number >
void LinearAlgebra::ReadWriteVector< Number >::import_elements ( const ::Vector< Number > &  vec,
VectorOperation::values  operation,
const std::shared_ptr< const Utilities::MPI::CommunicationPatternBase > &  communication_pattern = {} 
)

Imports all the elements present in the vector's IndexSet from the input vector vec. VectorOperation::values operation is used to decide if the elements in V should be added to the current vector or replace the current elements.

Note
The parameter communication_pattern is ignored since we are dealing with a serial vector here.

◆ import() [1/7]

template<typename Number >
void LinearAlgebra::ReadWriteVector< Number >::import ( const ::Vector< Number > &  V,
VectorOperation::values  operation,
const std::shared_ptr< const Utilities::MPI::CommunicationPatternBase > &  communication_pattern = {} 
)
inline
Deprecated:
Use import_elements() instead.

Definition at line 312 of file read_write_vector.h.

◆ import_elements() [2/9]

template<typename Number >
template<typename MemorySpace >
void LinearAlgebra::ReadWriteVector< Number >::import_elements ( const distributed::Vector< Number, MemorySpace > &  vec,
VectorOperation::values  operation,
const std::shared_ptr< const Utilities::MPI::CommunicationPatternBase > &  communication_pattern = {} 
)

Imports all the elements present in the vector's IndexSet from the input vector vec. VectorOperation::values operation is used to decide if the elements in V should be added to the current vector or replace the current elements. The last parameter can be used if the same communication pattern is used multiple times. This can be used to improve performance.

◆ import() [2/7]

template<typename Number >
template<typename MemorySpace >
void LinearAlgebra::ReadWriteVector< Number >::import ( const distributed::Vector< Number, MemorySpace > &  V,
VectorOperation::values  operation,
const std::shared_ptr< const Utilities::MPI::CommunicationPatternBase > &  communication_pattern = {} 
)
inline
Deprecated:
Use import_elements() instead.

Definition at line 341 of file read_write_vector.h.

◆ import_elements() [3/9]

template<typename Number >
void LinearAlgebra::ReadWriteVector< Number >::import_elements ( const PETScWrappers::MPI::Vector petsc_vec,
VectorOperation::values  operation,
const std::shared_ptr< const Utilities::MPI::CommunicationPatternBase > &  communication_pattern = {} 
)

Imports all the elements present in the vector's IndexSet from the input vector petsc_vec. VectorOperation::values operation is used to decide if the elements in V should be added to the current vector or replace the current elements. The last parameter can be used if the same communication pattern is used multiple times. This can be used to improve performance.

◆ import() [3/7]

template<typename Number >
void LinearAlgebra::ReadWriteVector< Number >::import ( const PETScWrappers::MPI::Vector V,
VectorOperation::values  operation,
const std::shared_ptr< const Utilities::MPI::CommunicationPatternBase > &  communication_pattern = {} 
)
inline
Deprecated:
Use import_elements() instead.

Definition at line 371 of file read_write_vector.h.

◆ import_elements() [4/9]

template<typename Number >
void LinearAlgebra::ReadWriteVector< Number >::import_elements ( const TrilinosWrappers::MPI::Vector trilinos_vec,
VectorOperation::values  operation,
const std::shared_ptr< const Utilities::MPI::CommunicationPatternBase > &  communication_pattern = {} 
)

Imports all the elements present in the vector's IndexSet from the input vector trilinos_vec. VectorOperation::values operation is used to decide if the elements in V should be added to the current vector or replace the current elements. The last parameter can be used if the same communication pattern is used multiple times. This can be used to improve performance.

Note
The trilinos_vec is not allowed to have ghost entries.

◆ import() [4/7]

template<typename Number >
void LinearAlgebra::ReadWriteVector< Number >::import ( const TrilinosWrappers::MPI::Vector V,
VectorOperation::values  operation,
const std::shared_ptr< const Utilities::MPI::CommunicationPatternBase > &  communication_pattern = {} 
)
inline
Deprecated:
Use import_elements() instead.

Definition at line 403 of file read_write_vector.h.

◆ import_elements() [5/9]

template<typename Number >
template<typename Dummy = Number>
std::enable_if_t< std::is_same_v< Dummy, Number > && ::is_tpetra_type< Number >::value > LinearAlgebra::ReadWriteVector< Number >::import_elements ( const TpetraWrappers::Vector< Number > &  tpetra_vec,
VectorOperation::values  operation,
const std::shared_ptr< const Utilities::MPI::CommunicationPatternBase > &  communication_pattern = {} 
)

Imports all the elements present in the vector's IndexSet from the input vector tpetra_vec. VectorOperation::values operation is used to decide if the elements in V should be added to the current vector or replace the current elements. The last parameter can be used if the same communication pattern is used multiple times. This can be used to improve performance.

◆ import() [5/7]

template<typename Number >
template<typename Dummy = Number>
std::enable_if_t< std::is_same_v< Dummy, Number > && ::is_tpetra_type< Number >::value > LinearAlgebra::ReadWriteVector< Number >::import ( const TpetraWrappers::Vector< Number > &  V,
VectorOperation::values  operation,
const std::shared_ptr< const Utilities::MPI::CommunicationPatternBase > &  communication_pattern = {} 
)
inline
Deprecated:
Use import_elements() instead.

Definition at line 435 of file read_write_vector.h.

◆ import_elements() [6/9]

template<typename Number >
void LinearAlgebra::ReadWriteVector< Number >::import_elements ( const EpetraWrappers::Vector epetra_vec,
VectorOperation::values  operation,
const std::shared_ptr< const Utilities::MPI::CommunicationPatternBase > &  communication_pattern = {} 
)

Imports all the elements present in the vector's IndexSet from the input vector epetra_vec. VectorOperation::values operation is used to decide if the elements in V should be added to the current vector or replace the current elements. The last parameter can be used if the same communication pattern is used multiple times. This can be used to improve performance.

◆ import() [6/7]

template<typename Number >
void LinearAlgebra::ReadWriteVector< Number >::import ( const EpetraWrappers::Vector V,
VectorOperation::values  operation,
const std::shared_ptr< const Utilities::MPI::CommunicationPatternBase > &  communication_pattern = {} 
)
inline
Deprecated:
Use import_elements() instead.

Definition at line 464 of file read_write_vector.h.

◆ import_elements() [7/9]

template<typename Number >
void LinearAlgebra::ReadWriteVector< Number >::import_elements ( const CUDAWrappers::Vector< Number > &  cuda_vec,
VectorOperation::values  operation,
const std::shared_ptr< const Utilities::MPI::CommunicationPatternBase > &  communication_pattern = {} 
)

Import all the elements present in the vector's IndexSet from the input vector cuda_vec. VectorOperation::values operation is used to decide if the elements in V should be added to the current vector or replace the current elements. The last parameter is not used.

◆ import() [7/7]

template<typename Number >
void LinearAlgebra::ReadWriteVector< Number >::import ( const CUDAWrappers::Vector< Number > &  V,
VectorOperation::values  operation,
const std::shared_ptr< const Utilities::MPI::CommunicationPatternBase > &  communication_pattern = {} 
)
inline
Deprecated:
Use import_elements() instead.

Definition at line 493 of file read_write_vector.h.

◆ size()

template<typename Number >
size_type LinearAlgebra::ReadWriteVector< Number >::size ( ) const
overridevirtual

The value returned by this function denotes the dimension of the vector spaces that are modeled by objects of this kind. However, objects of the current class do not actually stores all elements of vectors of this space but may, in fact store only a subset. The number of elements stored is returned by n_elements() and is smaller or equal to the number returned by the current function.

Implements ReadVector< Number >.

◆ locally_owned_size()

template<typename Number >
size_type LinearAlgebra::ReadWriteVector< Number >::locally_owned_size ( ) const

Return the local size of the vector, i.e., the number of indices owned locally.

◆ get_stored_elements()

template<typename Number >
const IndexSet & LinearAlgebra::ReadWriteVector< Number >::get_stored_elements ( ) const

Return the IndexSet that represents the indices of the elements stored.

◆ begin() [1/2]

template<typename Number >
iterator LinearAlgebra::ReadWriteVector< Number >::begin ( )

Make the ReadWriteVector class a bit like the vector<> class of the C++ standard library by returning iterators to the start and end of the locally stored elements of this vector.

◆ begin() [2/2]

template<typename Number >
const_iterator LinearAlgebra::ReadWriteVector< Number >::begin ( ) const

Return constant iterator to the start of the locally stored elements of the vector.

◆ end() [1/2]

template<typename Number >
iterator LinearAlgebra::ReadWriteVector< Number >::end ( )

Return an iterator pointing to the element past the end of the array of locally stored entries.

◆ end() [2/2]

template<typename Number >
const_iterator LinearAlgebra::ReadWriteVector< Number >::end ( ) const

Return a constant iterator pointing to the element past the end of the array of the locally stored entries.

◆ operator()() [1/2]

template<typename Number >
Number LinearAlgebra::ReadWriteVector< Number >::operator() ( const size_type  global_index) const

Read access to the data in the position corresponding to global_index. An exception is thrown if global_index is not stored by the current object.

◆ operator()() [2/2]

template<typename Number >
Number & LinearAlgebra::ReadWriteVector< Number >::operator() ( const size_type  global_index)

Read and write access to the data in the position corresponding to global_index. An exception is thrown if global_index is not stored by the current object.

◆ operator[]() [1/2]

template<typename Number >
Number LinearAlgebra::ReadWriteVector< Number >::operator[] ( const size_type  global_index) const

Read access to the data in the position corresponding to global_index. An exception is thrown if global_index is not stored by the current object.

This function does the same thing as operator().

◆ operator[]() [2/2]

template<typename Number >
Number & LinearAlgebra::ReadWriteVector< Number >::operator[] ( const size_type  global_index)

Read and write access to the data in the position corresponding to global_index. An exception is thrown if global_index is not stored by the current object.

This function does the same thing as operator().

◆ extract_subvector_to() [1/3]

template<typename Number >
template<typename Number2 >
void LinearAlgebra::ReadWriteVector< Number >::extract_subvector_to ( const std::vector< size_type > &  indices,
std::vector< Number2 > &  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

for (unsigned int i=0; i<indices.size(); ++i)
values[i] = v[indices[i]];
AlignedVector< Number > values
Precondition
The sizes of the indices and values arrays must be identical.

◆ extract_subvector_to() [2/3]

template<typename Number >
virtual void LinearAlgebra::ReadWriteVector< Number >::extract_subvector_to ( const ArrayView< const types::global_dof_index > &  indices,
ArrayView< Number > &  entries 
) const
overridevirtual

Extract a range of elements all at once.

Implements ReadVector< Number >.

◆ extract_subvector_to() [3/3]

template<typename Number >
template<typename ForwardIterator , typename OutputIterator >
void LinearAlgebra::ReadWriteVector< Number >::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

ForwardIterator indices_p = indices_begin;
OutputIterator values_p = values_begin;
while (indices_p != indices_end)
{
*values_p = v[*indices_p];
++indices_p;
++values_p;
}
Precondition
It must be possible to write into as many memory locations starting at values_begin as there are iterators between indices_begin and indices_end.

◆ local_element() [1/2]

template<typename Number >
Number LinearAlgebra::ReadWriteVector< Number >::local_element ( const size_type  local_index) const

Read access to the data field specified by local_index. When you access elements in the order in which they are stored, it is necessary that you know in which they are stored. In other words, you need to know the map between the global indices of the elements this class stores, and the local indices into the contiguous array of these global elements. For this, see the general documentation of this class.

Performance: Direct array access (fast).

◆ local_element() [2/2]

template<typename Number >
Number & LinearAlgebra::ReadWriteVector< Number >::local_element ( const size_type  local_index)

Read and write access to the data field specified by local_index. When you access elements in the order in which they are stored, it is necessary that you know in which they are stored. In other words, you need to know the map between the global indices of the elements this class stores, and the local indices into the contiguous array of these global elements. For this, see the general documentation of this class.

Performance: Direct array access (fast).

◆ add() [1/3]

template<typename Number >
template<typename Number2 >
void LinearAlgebra::ReadWriteVector< Number >::add ( const std::vector< size_type > &  indices,
const std::vector< Number2 > &  values 
)

This function adds a whole set of values stored in values to the vector components specified by indices.

◆ add() [2/3]

template<typename Number >
template<typename Number2 >
void LinearAlgebra::ReadWriteVector< Number >::add ( const std::vector< size_type > &  indices,
const ReadWriteVector< Number2 > &  values 
)

This function is similar to the previous one but takes a ReadWriteVector of values.

◆ add() [3/3]

template<typename Number >
template<typename Number2 >
void LinearAlgebra::ReadWriteVector< Number >::add ( const size_type  n_elements,
const size_type indices,
const Number2 *  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.

◆ print()

template<typename Number >
void LinearAlgebra::ReadWriteVector< Number >::print ( std::ostream &  out,
const unsigned int  precision = 3,
const bool  scientific = true 
) const

Prints the vector to the output stream out.

◆ memory_consumption()

template<typename Number >
std::size_t LinearAlgebra::ReadWriteVector< Number >::memory_consumption ( ) const

Return the memory consumption of this class in bytes.

◆ import_elements() [8/9]

template<typename Number >
template<typename Dummy = Number>
std::enable_if_t< std::is_same_v< Dummy, Number > && ::is_tpetra_type< Number >::value > LinearAlgebra::ReadWriteVector< Number >::import_elements ( const Tpetra::Vector< Number, int, types::signed_global_dof_index > &  tpetra_vector,
const IndexSet locally_owned_elements,
VectorOperation::values  operation,
const MPI_Comm  mpi_comm,
const std::shared_ptr< const Utilities::MPI::CommunicationPatternBase > &  communication_pattern 
)
protected

Import all the elements present in the vector's IndexSet from the input vector tpetra_vector. This is an helper function and it should not be used directly.

◆ import_elements() [9/9]

template<typename Number >
void LinearAlgebra::ReadWriteVector< Number >::import_elements ( const Epetra_MultiVector &  multivector,
const IndexSet locally_owned_elements,
VectorOperation::values  operation,
const MPI_Comm  mpi_comm,
const std::shared_ptr< const Utilities::MPI::CommunicationPatternBase > &  communication_pattern 
)
protected

Import all the elements present in the vector's IndexSet from the input vector multivector. This is an helper function and it should not be used directly.

◆ global_to_local()

template<typename Number >
unsigned int LinearAlgebra::ReadWriteVector< Number >::global_to_local ( const types::global_dof_index  global_index) const
protected

Return the local position of global_index.

◆ create_tpetra_comm_pattern()

template<typename Number >
TpetraWrappers::CommunicationPattern LinearAlgebra::ReadWriteVector< Number >::create_tpetra_comm_pattern ( const IndexSet source_index_set,
const MPI_Comm  mpi_comm 
)
protected

Return a TpetraWrappers::CommunicationPattern and store it for future use.

◆ create_epetra_comm_pattern()

template<typename Number >
EpetraWrappers::CommunicationPattern LinearAlgebra::ReadWriteVector< Number >::create_epetra_comm_pattern ( const IndexSet source_index_set,
const MPI_Comm  mpi_comm 
)
protected

Return a EpetraWrappers::CommunicationPattern and store it for future use.

◆ subscribe()

void Subscriptor::subscribe ( std::atomic< bool > *const  validity,
const std::string &  identifier = "" 
) const
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 135 of file subscriptor.cc.

◆ unsubscribe()

void Subscriptor::unsubscribe ( std::atomic< bool > *const  validity,
const std::string &  identifier = "" 
) const
inherited

Unsubscribes a user from the object.

Note
The identifier and the validity pointer must be the same as the one supplied to subscribe().

Definition at line 155 of file subscriptor.cc.

◆ n_subscriptions()

unsigned int Subscriptor::n_subscriptions ( ) const
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 300 of file subscriptor.h.

◆ list_subscribers() [1/2]

template<typename StreamType >
void Subscriptor::list_subscribers ( StreamType &  stream) const
inlineinherited

List the subscribers to the input stream.

Definition at line 317 of file subscriptor.h.

◆ list_subscribers() [2/2]

void Subscriptor::list_subscribers ( ) const
inherited

List the subscribers to deallog.

Definition at line 203 of file subscriptor.cc.

◆ serialize()

template<class Archive >
void Subscriptor::serialize ( Archive &  ar,
const unsigned int  version 
)
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 309 of file subscriptor.h.

◆ check_no_subscribers()

void Subscriptor::check_no_subscribers ( ) const
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.

Note
Since this function is just a consistency check it does nothing in release mode.
If this function is called when there is an uncaught exception then, rather than aborting, this function prints an error message to the standard error stream and returns.

Definition at line 52 of file subscriptor.cc.

Friends And Related Symbol Documentation

◆ ReadWriteVector

template<typename Number >
template<typename Number2 >
friend class ReadWriteVector
friend

Definition at line 828 of file read_write_vector.h.

Member Data Documentation

◆ stored_elements

template<typename Number >
IndexSet LinearAlgebra::ReadWriteVector< Number >::stored_elements
protected

Indices of the elements stored.

Definition at line 801 of file read_write_vector.h.

◆ source_stored_elements

template<typename Number >
IndexSet LinearAlgebra::ReadWriteVector< Number >::source_stored_elements
protected

IndexSet of the elements of the last imported vector;

Definition at line 806 of file read_write_vector.h.

◆ comm_pattern

template<typename Number >
std::shared_ptr<Utilities::MPI::CommunicationPatternBase> LinearAlgebra::ReadWriteVector< Number >::comm_pattern
protected

CommunicationPattern for the communication between the source_stored_elements IndexSet and the current vector.

Definition at line 812 of file read_write_vector.h.

◆ values

template<typename Number >
AlignedVector<Number> LinearAlgebra::ReadWriteVector< Number >::values
protected

Locally stored elements.

Definition at line 817 of file read_write_vector.h.

◆ thread_loop_partitioner

template<typename Number >
std::shared_ptr<::parallel::internal::TBBPartitioner> LinearAlgebra::ReadWriteVector< Number >::thread_loop_partitioner
mutableprotected

For parallel loops with TBB, this member variable stores the affinity information of loops.

Definition at line 824 of file read_write_vector.h.

◆ counter

std::atomic<unsigned int> Subscriptor::counter
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 218 of file subscriptor.h.

◆ counter_map

std::map<std::string, unsigned int> Subscriptor::counter_map
mutableprivateinherited

In this map, we count subscriptions for each different identification string supplied to subscribe().

Definition at line 224 of file subscriptor.h.

◆ validity_pointers

std::vector<std::atomic<bool> *> Subscriptor::validity_pointers
mutableprivateinherited

In this vector, we store pointers to the validity bool in the SmartPointer objects that subscribe to this class.

Definition at line 240 of file subscriptor.h.

◆ object_info

const std::type_info* Subscriptor::object_info
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 248 of file subscriptor.h.

◆ mutex

std::mutex Subscriptor::mutex
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 271 of file subscriptor.h.


The documentation for this class was generated from the following file: