Reference documentation for deal.II version 9.5.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\}}\)
Loading...
Searching...
No Matches
Public Types | Public Member Functions | Static Public Member Functions | Private Types | Private Member Functions | Private Attributes | Static Private Attributes | Friends | Related Symbols | List of all members
Vector< Number > Class Template Reference

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

Inheritance diagram for Vector< Number >:
[legend]

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

Vector< int >::real_type lp_norm (const real_type) const
 
Basic object handling
 Vector ()
 
 Vector (const Vector< Number > &v)
 
 Vector (Vector< Number > &&v) noexcept=default
 
template<typename OtherNumber >
 Vector (const Vector< OtherNumber > &v)
 
template<typename OtherNumber >
 Vector (const std::initializer_list< OtherNumber > &v)
 
 Vector (const PETScWrappers::VectorBase &v)
 
 Vector (const TrilinosWrappers::MPI::Vector &v)
 
 Vector (const size_type n)
 
template<typename InputIterator >
 Vector (const InputIterator first, const InputIterator last)
 
virtual ~Vector () override=default
 
void compress (VectorOperation::values operation=VectorOperation::unknown) const
 
virtual void reinit (const size_type N, const bool omit_zeroing_entries=false)
 
void grow_or_shrink (const size_type N)
 
void apply_givens_rotation (const std::array< Number, 3 > &csr, const size_type i, const size_type k)
 
template<typename Number2 >
void reinit (const Vector< Number2 > &V, const bool omit_zeroing_entries=false)
 
virtual void swap (Vector< Number > &v)
 
Vector< Number > & operator= (const Number s)
 
Vector< Number > & operator= (const Vector< Number > &v)
 
Vector< Number > & operator= (Vector< Number > &&v) noexcept=default
 
template<typename Number2 >
Vector< Number > & operator= (const Vector< Number2 > &v)
 
Vector< Number > & operator= (const BlockVector< Number > &v)
 
Vector< Number > & operator= (const PETScWrappers::VectorBase &v)
 
Vector< Number > & operator= (const TrilinosWrappers::MPI::Vector &v)
 
template<typename Number2 >
bool operator== (const Vector< Number2 > &v) const
 
template<typename Number2 >
bool operator!= (const Vector< Number2 > &v) const
 
Scalar products, norms and related operations
template<typename Number2 >
Number operator* (const Vector< Number2 > &V) const
 
real_type norm_sqr () const
 
Number mean_value () const
 
real_type l1_norm () const
 
real_type l2_norm () const
 
real_type lp_norm (const real_type p) const
 
real_type linfty_norm () const
 
Number add_and_dot (const Number a, const Vector< Number > &V, const Vector< Number > &W)
 
Data access
pointer data ()
 
const_pointer data () const
 
iterator begin ()
 
const_iterator begin () const
 
iterator end ()
 
const_iterator end () const
 
Number operator() (const size_type i) const
 
Number & operator() (const size_type i)
 
Number operator[] (const size_type i) const
 
Number & operator[] (const size_type i)
 
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
 
Modification of vectors
Vector< Number > & operator+= (const Vector< Number > &V)
 
Vector< Number > & operator-= (const Vector< Number > &V)
 
template<typename OtherNumber >
void add (const std::vector< size_type > &indices, const std::vector< OtherNumber > &values)
 
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 add (const Number s)
 
void add (const Number a, const Vector< Number > &V, const Number b, const Vector< Number > &W)
 
void add (const Number a, const Vector< Number > &V)
 
void sadd (const Number s, const Vector< Number > &V)
 
void sadd (const Number s, const Number a, const Vector< Number > &V)
 
Vector< Number > & operator*= (const Number factor)
 
Vector< Number > & operator/= (const Number factor)
 
void scale (const Vector< Number > &scaling_factors)
 
template<typename Number2 >
void scale (const Vector< Number2 > &scaling_factors)
 
void equ (const Number a, const Vector< Number > &u)
 
template<typename Number2 >
void equ (const Number a, const Vector< Number2 > &u)
 
void update_ghost_values () const
 
Input and output
void print (std::ostream &out, const unsigned int precision=3, const bool scientific=true, const bool across=true) const
 
void block_write (std::ostream &out) const
 
void block_read (std::istream &in)
 
template<class Archive >
void save (Archive &ar, const unsigned int version) const
 
template<class Archive >
void load (Archive &ar, const unsigned int version)
 
template<class Archive >
void serialize (Archive &archive, const unsigned int version)
 
Information about the object
bool in_local_range (const size_type global_index) const
 
IndexSet locally_owned_elements () const
 
size_type size () const
 
size_type locally_owned_size () const
 
bool all_zero () const
 
bool is_non_negative () const
 
std::size_t memory_consumption () const
 
bool has_ghost_elements () const
 
void zero_out_ghost_values () 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)
 

Private Types

using map_value_type = decltype(counter_map)::value_type
 
using map_iterator = decltype(counter_map)::iterator
 

Private Member Functions

void maybe_reset_thread_partitioner ()
 
void do_reinit (const size_type new_size, const bool omit_zeroing_entries, const bool reset_partitioner)
 
void check_no_subscribers () const noexcept
 

Private Attributes

AlignedVector< Number > values
 
std::shared_ptr< parallel::internal::TBBPartitionerthread_loop_partitioner
 
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 Vector
 

Related Symbols

(Note that these are not member symbols.)

template<typename Number >
void swap (LinearAlgebra::CUDAWrappers::Vector< Number > &u, LinearAlgebra::CUDAWrappers::Vector< Number > &v)
 
template<typename Number , typename MemorySpace >
void swap (LinearAlgebra::distributed::Vector< Number, MemorySpace > &u, LinearAlgebra::distributed::Vector< Number, MemorySpace > &v)
 
template<typename Number >
void swap (LinearAlgebra::ReadWriteVector< Number > &u, LinearAlgebra::ReadWriteVector< Number > &v)
 
template<typename Number >
void swap (Vector< Number > &u, Vector< Number > &v)
 
template<typename number >
std::ostream & operator<< (std::ostream &out, const Vector< number > &v)
 

Detailed Description

template<typename Number>
class Vector< Number >

A class that represents a vector of numerical elements. As for the other classes, in the Vector classes group, this class has a substantial number of member functions. These include:

In contrast to the C++ standard library class std::vector, this class intends to implement not simply an array that allows access to its elements, but indeed a vector that is a member of the mathematical concept of a "vector space" suitable for numerical computations.

Note
Instantiations for this template are provided for <float>, <double>, <std::complex<float>>, <std::complex<double>>; others can be generated in application programs (see the section on Template instantiations in the manual).

Definition at line 108 of file vector.h.

Member Typedef Documentation

◆ value_type

template<typename Number >
using Vector< Number >::value_type = Number

This class only supports basic numeric types (i.e., we support double and float but not automatically differentiated numbers).

Note
we test real_type here to get the underlying scalar type when using std::complex. Declare standard types used in all containers. These types parallel those in the C++ standard libraries vector<...> class.

Definition at line 128 of file vector.h.

◆ pointer

template<typename Number >
using Vector< Number >::pointer = value_type *

Definition at line 129 of file vector.h.

◆ const_pointer

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

Definition at line 130 of file vector.h.

◆ iterator

template<typename Number >
using Vector< Number >::iterator = value_type *

Definition at line 131 of file vector.h.

◆ const_iterator

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

Definition at line 132 of file vector.h.

◆ reference

template<typename Number >
using Vector< Number >::reference = value_type &

Definition at line 133 of file vector.h.

◆ const_reference

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

Definition at line 134 of file vector.h.

◆ size_type

template<typename Number >
using Vector< Number >::size_type = types::global_dof_index

Definition at line 135 of file vector.h.

◆ real_type

template<typename Number >
using Vector< Number >::real_type = typename numbers::NumberTraits<Number>::real_type

Declare a type that has holds real-valued numbers with the same precision as the template argument to this class. If the template argument of this class is a real data type, then real_type equals the template argument. If the template argument is a std::complex type then real_type equals the type underlying the complex numbers.

This alias is used to represent the return type of norms.

Definition at line 146 of file 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 230 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 235 of file subscriptor.h.

Constructor & Destructor Documentation

◆ Vector() [1/9]

template<typename Number >
Vector< Number >::Vector ( )

Constructor. Create a vector of dimension zero.

◆ Vector() [2/9]

template<typename Number >
Vector< Number >::Vector ( const Vector< Number > &  v)

Copy constructor. Sets the dimension to that of the given vector, and copies all elements.

We would like to make this constructor explicit, but standard containers insist on using it implicitly.

Note
If deal.II is configured with threads, this operation will run multi-threaded by splitting the work into smaller chunks (assuming there is enough work to make this worthwhile).

◆ Vector() [3/9]

template<typename Number >
Vector< Number >::Vector ( Vector< Number > &&  v)
defaultnoexcept

Move constructor. Creates a new vector by stealing the internal data of the vector v.

◆ Vector() [4/9]

template<typename Number >
template<typename OtherNumber >
Vector< Number >::Vector ( const Vector< OtherNumber > &  v)
explicit

Copy constructor taking a vector of another data type.

This constructor will fail to compile if there is no conversion path from OtherNumber to Number. You may lose accuracy when copying to a vector with data elements with less accuracy.

◆ Vector() [5/9]

template<typename Number >
template<typename OtherNumber >
Vector< Number >::Vector ( const std::initializer_list< OtherNumber > &  v)
explicit

Copy constructor taking an object of type std::initializer_list. This constructor can be used to initialize a vector using a brace-enclosed list of numbers, such as in the following example:

Vector<double> v({1,2,3});

This creates a vector of size 3, whose (double precision) elements have values 1.0, 2.0, and 3.0.

This constructor will fail to compile if there is no conversion path from OtherNumber to Number. You may lose accuracy when copying to a vector with data elements with less accuracy.

◆ Vector() [6/9]

template<typename Number >
Vector< Number >::Vector ( const PETScWrappers::VectorBase v)
explicit

Another copy constructor: copy the values from a PETSc vector class. This copy constructor is only available if PETSc was detected during configuration time.

Note that due to the communication model used in MPI, this operation can only succeed if all processes do it at the same time when v is a distributed vector: It is not possible for only one process to obtain a copy of a parallel vector while the other jobs do something else.

◆ Vector() [7/9]

template<typename Number >
Vector< Number >::Vector ( const TrilinosWrappers::MPI::Vector< Number > &  v)
explicit

Another copy constructor: copy the values from a Trilinos wrapper vector. This copy constructor is only available if Trilinos was detected during configuration time.

Note
Due to the communication model used in MPI, this operation can only succeed if all processes that have knowledge of v (i.e. those given by v.get_mpi_communicator()) do it at the same time. This means that unless you use a split MPI communicator then it is not normally possible for only one or a subset of processes to obtain a copy of a parallel vector while the other jobs do something else. In other words, calling this function is a 'collective operation' that needs to be executed by all MPI processes that jointly share v.

◆ Vector() [8/9]

template<typename Number >
Vector< Number >::Vector ( const size_type  n)
explicit

Constructor. Set dimension to n and initialize all elements with zero.

The constructor is made explicit to avoid accidents like this: v=0;. Presumably, the user wants to set every element of the vector to zero, but instead, what happens is this call: v=Vector<number>(0);, i.e. the vector is replaced by one of length zero.

◆ Vector() [9/9]

template<typename Number >
template<typename InputIterator >
Vector< Number >::Vector ( const InputIterator  first,
const InputIterator  last 
)

Initialize the vector with a given range of values pointed to by the iterators. This function is there in analogy to the std::vector class.

◆ ~Vector()

template<typename Number >
virtual Vector< Number >::~Vector ( )
overridevirtualdefault

Destructor, deallocates memory. Made virtual to allow for derived classes to behave properly.

Member Function Documentation

◆ compress()

template<typename Number >
void Vector< Number >::compress ( VectorOperation::values  operation = VectorOperation::unknown) const

This function does nothing but exists for compatibility with the parallel vector classes.

For the parallel vector wrapper class, this function compresses the underlying representation of the vector, i.e. flushes the buffers of the vector object if it has any. This function is necessary after writing into a vector element-by-element and before anything else can be done on it.

However, for the implementation of this class, it is immaterial and thus an empty function.

◆ reinit() [1/2]

template<typename Number >
virtual void Vector< Number >::reinit ( const size_type  N,
const bool  omit_zeroing_entries = false 
)
virtual

Change the dimension of the vector to N. The reserved memory for this vector remains unchanged if possible, to make things faster; this may waste some memory, so keep this in mind. However, if N==0 all memory is freed, i.e. if you want to resize the vector and release the memory not needed, you have to first call reinit(0) and then reinit(N). This cited behavior is analogous to that of the standard library containers.

If omit_zeroing_entries is false, the vector is filled by zeros. Otherwise, the elements are left an unspecified state.

This function is virtual in order to allow for derived classes to handle memory separately.

◆ grow_or_shrink()

template<typename Number >
void Vector< Number >::grow_or_shrink ( const size_type  N)

Same as above, but will preserve the values of vector upon resizing. If we new size is bigger, we will have

\[ \mathbf V \rightarrow \left( \begin{array}{c} \mathbf V \\ \mathbf 0 \end{array} \right) \]

whereas if the desired size is smaller, then

\[ \left( \begin{array}{c} \mathbf V_1 \\ \mathbf V_2 \end{array} \right) \rightarrow \mathbf V_1 \]

◆ apply_givens_rotation()

template<typename Number >
void Vector< Number >::apply_givens_rotation ( const std::array< Number, 3 > &  csr,
const size_type  i,
const size_type  k 
)

Apply Givens rotation csr (a triplet of cosine, sine and radius, see Utilities::LinearAlgebra::givens_rotation()) to the vector in the plane spanned by the i'th and k'th unit vectors.

◆ reinit() [2/2]

template<typename Number >
template<typename Number2 >
void Vector< Number >::reinit ( const Vector< Number2 > &  V,
const bool  omit_zeroing_entries = false 
)

Change the dimension to that of the vector V. The same applies as for the other reinit function.

The elements of V are not copied, i.e. this function is the same as calling reinit (V.size(), omit_zeroing_entries).

◆ swap()

template<typename Number >
virtual void Vector< Number >::swap ( Vector< Number > &  v)
virtual

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.

This function is virtual in order to allow for derived classes to handle memory separately.

◆ operator=() [1/7]

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

Set all components of the vector to the given number s.

Since the semantics of assigning a scalar to a vector are not immediately clear, this operator should really only be used if you want to set the entire vector to zero. This allows the intuitive notation v=0. Assigning other values is deprecated and may be disallowed in the future.

Note
If deal.II is configured with threads, this operation will run multi-threaded by splitting the work into smaller chunks (assuming there is enough work to make this worthwhile).

◆ operator=() [2/7]

template<typename Number >
Vector< Number > & Vector< Number >::operator= ( const Vector< Number > &  v)

Copy the given vector. Resize the present vector if necessary.

Note
If deal.II is configured with threads, this operation will run multi-threaded by splitting the work into smaller chunks (assuming there is enough work to make this worthwhile).

◆ operator=() [3/7]

template<typename Number >
Vector< Number > & Vector< Number >::operator= ( Vector< Number > &&  v)
defaultnoexcept

Move the given vector. This operator replaces the present vector with the internal data of the vector v and resets v to the state it would have after being newly default-constructed.

◆ operator=() [4/7]

template<typename Number >
template<typename Number2 >
Vector< Number > & Vector< Number >::operator= ( const Vector< Number2 > &  v)

Copy the given vector. Resize the present vector if necessary.

Note
If deal.II is configured with threads, this operation will run multi-threaded by splitting the work into smaller chunks (assuming there is enough work to make this worthwhile).

◆ operator=() [5/7]

template<typename Number >
Vector< Number > & Vector< Number >::operator= ( const BlockVector< Number > &  v)

Copy operator for assigning a block vector to a regular vector.

◆ operator=() [6/7]

template<typename Number >
Vector< Number > & Vector< Number >::operator= ( const PETScWrappers::VectorBase v)

Another copy operator: copy the values from a PETSc wrapper vector class. This operator is only available if PETSc was detected during configuration time.

Note that due to the communication model used in MPI, this operation can only succeed if all processes do it at the same time when v is a distributed vector: It is not possible for only one process to obtain a copy of a parallel vector while the other jobs do something else.

◆ operator=() [7/7]

template<typename Number >
Vector< Number > & Vector< Number >::operator= ( const TrilinosWrappers::MPI::Vector< Number > &  v)

Another copy operator: copy the values from a (sequential or parallel, depending on the underlying compiler) Trilinos wrapper vector class. This operator is only available if Trilinos was detected during configuration time.

Note
Due to the communication model used in MPI, this operation can only succeed if all processes that have knowledge of v (i.e. those given by v.get_mpi_communicator()) do it at the same time. This means that unless you use a split MPI communicator then it is not normally possible for only one or a subset of processes to obtain a copy of a parallel vector while the other jobs do something else. In other words, calling this function is a 'collective operation' that needs to be executed by all MPI processes that jointly share v.

◆ operator==()

template<typename Number >
template<typename Number2 >
bool Vector< Number >::operator== ( const Vector< Number2 > &  v) const

Test for equality. This function assumes that the present vector and the one to compare with have the same size already, since comparing vectors of different sizes makes not much sense anyway.

◆ operator!=()

template<typename Number >
template<typename Number2 >
bool Vector< Number >::operator!= ( const Vector< Number2 > &  v) const

Test for inequality. This function assumes that the present vector and the one to compare with have the same size already, since comparing vectors of different sizes makes not much sense anyway.

◆ operator*()

template<typename Number >
template<typename Number2 >
Number Vector< Number >::operator* ( const Vector< Number2 > &  V) const

Return the scalar product of two vectors. The return type is the underlying type of this vector, so the return type and the accuracy with which it the result is computed depend on the order of the arguments of this vector.

For complex vectors, the scalar product is implemented as \(\left<v,w\right>=\sum_i v_i \bar{w_i}\).

Note
If deal.II is configured with threads, this operation will run multi-threaded by splitting the work into smaller chunks (assuming there is enough work to make this worthwhile). The algorithm uses pairwise summation with the same order of summation in every run, which gives fully repeatable results from one run to another.

◆ norm_sqr()

template<typename Number >
real_type Vector< Number >::norm_sqr ( ) const

Return the square of the \(l_2\)-norm.

Note
If deal.II is configured with threads, this operation will run multi-threaded by splitting the work into smaller chunks (assuming there is enough work to make this worthwhile). The algorithm uses pairwise summation with the same order of summation in every run, which gives fully repeatable results from one run to another.

◆ mean_value()

template<typename Number >
Number Vector< Number >::mean_value ( ) const

Mean value of the elements of this vector.

Note
If deal.II is configured with threads, this operation will run multi-threaded by splitting the work into smaller chunks (assuming there is enough work to make this worthwhile). The algorithm uses pairwise summation with the same order of summation in every run, which gives fully repeatable results from one run to another.

◆ l1_norm()

template<typename Number >
real_type Vector< Number >::l1_norm ( ) const

\(l_1\)-norm of the vector. The sum of the absolute values.

Note
If deal.II is configured with threads, this operation will run multi-threaded by splitting the work into smaller chunks (assuming there is enough work to make this worthwhile). The algorithm uses pairwise summation with the same order of summation in every run, which gives fully repeatable results from one run to another.

◆ l2_norm()

template<typename Number >
real_type Vector< Number >::l2_norm ( ) const

\(l_2\)-norm of the vector. The square root of the sum of the squares of the elements.

Note
If deal.II is configured with threads, this operation will run multi-threaded by splitting the work into smaller chunks (assuming there is enough work to make this worthwhile). The algorithm uses pairwise summation with the same order of summation in every run, which gives fully repeatable results from one run to another.

◆ lp_norm() [1/2]

template<typename Number >
real_type Vector< Number >::lp_norm ( const real_type  p) const

\(l_p\)-norm of the vector. The pth root of the sum of the pth powers of the absolute values of the elements.

Note
If deal.II is configured with threads, this operation will run multi-threaded by splitting the work into smaller chunks (assuming there is enough work to make this worthwhile). The algorithm uses pairwise summation with the same order of summation in every run, which gives fully repeatable results from one run to another.

◆ linfty_norm()

template<typename Number >
real_type Vector< Number >::linfty_norm ( ) const

Maximum absolute value of the elements.

◆ add_and_dot()

template<typename Number >
Number Vector< Number >::add_and_dot ( const Number  a,
const Vector< Number > &  V,
const Vector< Number > &  W 
)

Performs a combined operation of a vector addition and a subsequent inner product, returning the value of the inner product. In other words, the result of this function is the same as if the user called

this->add(a, V);
return_value = *this * W;
void add(const std::vector< size_type > &indices, const std::vector< OtherNumber > &values)

The reason this function exists is that this operation involves less memory transfer than calling the two functions separately. This method only needs to load three vectors, this, V, W, whereas calling separate methods means to load the calling vector this twice. Since most vector operations are memory transfer limited, this reduces the time by 25% (or 50% if W equals this).

For complex-valued vectors, the scalar product in the second step is implemented as \(\left<v,w\right>=\sum_i v_i \bar{w_i}\).

Note
If deal.II is configured with threads, this operation will run multi-threaded by splitting the work into smaller chunks (assuming there is enough work to make this worthwhile). The algorithm uses pairwise summation with the same order of summation in every run, which gives fully repeatable results from one run to another.

◆ data() [1/2]

template<typename Number >
pointer Vector< Number >::data ( )

Return a pointer to the underlying data buffer.

◆ data() [2/2]

template<typename Number >
const_pointer Vector< Number >::data ( ) const

Return a const pointer to the underlying data buffer.

◆ begin() [1/2]

template<typename Number >
iterator Vector< Number >::begin ( )

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

◆ begin() [2/2]

template<typename Number >
const_iterator Vector< Number >::begin ( ) const

Return constant iterator to the start of the vectors.

◆ end() [1/2]

template<typename Number >
iterator Vector< Number >::end ( )

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

◆ end() [2/2]

template<typename Number >
const_iterator Vector< Number >::end ( ) const

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

◆ operator()() [1/2]

template<typename Number >
Number Vector< Number >::operator() ( const size_type  i) const

Access the value of the ith component.

◆ operator()() [2/2]

template<typename Number >
Number & Vector< Number >::operator() ( const size_type  i)

Access the ith component as a writeable reference.

◆ operator[]() [1/2]

template<typename Number >
Number Vector< Number >::operator[] ( const size_type  i) const

Access the value of the ith component.

Exactly the same as operator().

◆ operator[]() [2/2]

template<typename Number >
Number & Vector< Number >::operator[] ( const size_type  i)

Access the ith component as a writeable reference.

Exactly the same as operator().

◆ extract_subvector_to() [1/2]

template<typename Number >
template<typename OtherNumber >
void Vector< Number >::extract_subvector_to ( const std::vector< size_type > &  indices,
std::vector< OtherNumber > &  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
Definition vector.h:1022
Precondition
The sizes of the indices and values arrays must be identical.

◆ extract_subvector_to() [2/2]

template<typename Number >
template<typename ForwardIterator , typename OutputIterator >
void Vector< 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.

◆ operator+=()

template<typename Number >
Vector< Number > & Vector< Number >::operator+= ( const Vector< Number > &  V)

Add the given vector to the present one.

Note
If deal.II is configured with threads, this operation will run multi-threaded by splitting the work into smaller chunks (assuming there is enough work to make this worthwhile).

◆ operator-=()

template<typename Number >
Vector< Number > & Vector< Number >::operator-= ( const Vector< Number > &  V)

Subtract the given vector from the present one.

Note
If deal.II is configured with threads, this operation will run multi-threaded by splitting the work into smaller chunks (assuming there is enough work to make this worthwhile).

◆ add() [1/6]

template<typename Number >
template<typename OtherNumber >
void Vector< Number >::add ( const std::vector< size_type > &  indices,
const std::vector< OtherNumber > &  values 
)

A collective add operation: This function adds a whole set of values stored in values to the vector components specified by indices.

◆ add() [2/6]

template<typename Number >
template<typename OtherNumber >
void Vector< Number >::add ( const std::vector< size_type > &  indices,
const Vector< OtherNumber > &  values 
)

This is a second collective add operation. As a difference, this function takes a deal.II vector of values.

◆ add() [3/6]

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

◆ add() [4/6]

template<typename Number >
void Vector< Number >::add ( const Number  s)

Addition of s to all components. Note that s is a scalar and not a vector.

Note
If deal.II is configured with threads, this operation will run multi-threaded by splitting the work into smaller chunks (assuming there is enough work to make this worthwhile).

◆ add() [5/6]

template<typename Number >
void Vector< Number >::add ( const Number  a,
const Vector< Number > &  V,
const Number  b,
const Vector< Number > &  W 
)

Multiple addition of scaled vectors, i.e. *this += a*V+b*W.

Note
If deal.II is configured with threads, this operation will run multi-threaded by splitting the work into smaller chunks (assuming there is enough work to make this worthwhile).

◆ add() [6/6]

template<typename Number >
void Vector< Number >::add ( const Number  a,
const Vector< Number > &  V 
)

Simple addition of a multiple of a vector, i.e. *this += a*V.

Note
If deal.II is configured with threads, this operation will run multi-threaded by splitting the work into smaller chunks (assuming there is enough work to make this worthwhile).

◆ sadd() [1/2]

template<typename Number >
void Vector< Number >::sadd ( const Number  s,
const Vector< Number > &  V 
)

Scaling and simple vector addition, i.e. this = s(*this)+V.

Note
If deal.II is configured with threads, this operation will run multi-threaded by splitting the work into smaller chunks (assuming there is enough work to make this worthwhile).

◆ sadd() [2/2]

template<typename Number >
void Vector< Number >::sadd ( const Number  s,
const Number  a,
const Vector< Number > &  V 
)

Scaling and simple addition, i.e. this = s(*this)+a*V.

Note
If deal.II is configured with threads, this operation will run multi-threaded by splitting the work into smaller chunks (assuming there is enough work to make this worthwhile).

◆ operator*=()

template<typename Number >
Vector< Number > & Vector< Number >::operator*= ( const Number  factor)

Scale each element of the vector by a constant value.

Note
If deal.II is configured with threads, this operation will run multi-threaded by splitting the work into smaller chunks (assuming there is enough work to make this worthwhile).

◆ operator/=()

template<typename Number >
Vector< Number > & Vector< Number >::operator/= ( const Number  factor)

Scale each element of the vector by the inverse of the given value.

Note
If deal.II is configured with threads, this operation will run multi-threaded by splitting the work into smaller chunks (assuming there is enough work to make this worthwhile).

◆ scale() [1/2]

template<typename Number >
void Vector< Number >::scale ( const Vector< Number > &  scaling_factors)

Scale each element of this vector by the corresponding element in the argument. This function is mostly meant to simulate multiplication (and immediate re-assignment) by a diagonal scaling matrix.

Note
If deal.II is configured with threads, this operation will run multi-threaded by splitting the work into smaller chunks (assuming there is enough work to make this worthwhile).

◆ scale() [2/2]

template<typename Number >
template<typename Number2 >
void Vector< Number >::scale ( const Vector< Number2 > &  scaling_factors)

Scale each element of this vector by the corresponding element in the argument. This function is mostly meant to simulate multiplication (and immediate re-assignment) by a diagonal scaling matrix.

◆ equ() [1/2]

template<typename Number >
void Vector< Number >::equ ( const Number  a,
const Vector< Number > &  u 
)

Assignment *this = a*u.

Note
If deal.II is configured with threads, this operation will run multi-threaded by splitting the work into smaller chunks (assuming there is enough work to make this worthwhile).

◆ equ() [2/2]

template<typename Number >
template<typename Number2 >
void Vector< Number >::equ ( const Number  a,
const Vector< Number2 > &  u 
)

Assignment *this = a*u.

◆ update_ghost_values()

template<typename Number >
void Vector< Number >::update_ghost_values ( ) const

This function does nothing but exists for compatibility with the parallel vector classes (e.g., LinearAlgebra::distributed::Vector class).

◆ print()

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

Print to a stream. precision denotes the desired precision with which values shall be printed, scientific whether scientific notation shall be used. If across is true then the vector is printed in a line, while if false then the elements are printed on a separate line each.

◆ block_write()

template<typename Number >
void Vector< Number >::block_write ( std::ostream &  out) const

Write the vector en bloc to a file. This is done in a binary mode, so the output is neither readable by humans nor (probably) by other computers using a different operating system or number format.

◆ block_read()

template<typename Number >
void Vector< Number >::block_read ( std::istream &  in)

Read a vector en block from a file. This is done using the inverse operations to the above function, so it is reasonably fast because the bitstream is not interpreted.

The vector is resized if necessary.

A primitive form of error checking is performed which will recognize the bluntest attempts to interpret some data as a vector stored bitwise to a file, but not more.

◆ save()

template<typename Number >
template<class Archive >
void Vector< Number >::save ( Archive &  ar,
const unsigned int  version 
) const

Write the data of this object to a stream for the purpose of serialization using the BOOST serialization library.

◆ load()

template<typename Number >
template<class Archive >
void Vector< Number >::load ( Archive &  ar,
const unsigned int  version 
)

Read the data of this object from a stream for the purpose of serialization using the BOOST serialization library.

◆ serialize()

template<typename Number >
template<class Archive >
void Vector< Number >::serialize ( Archive &  archive,
const unsigned int  version 
)

Write and read the data of this object from a stream for the purpose of serialization using the BOOST serialization library.

◆ in_local_range()

template<typename Number >
bool Vector< Number >::in_local_range ( const size_type  global_index) const

Return true if the given global index is in the local range of this processor. Since this is not a distributed vector the method always returns true.

◆ locally_owned_elements()

template<typename Number >
IndexSet Vector< Number >::locally_owned_elements ( ) const

Return an index set that describes which elements of this vector are owned by the current processor. Note that this index set does not include elements this vector may store locally as ghost elements but that are in fact owned by another processor. As a consequence, the index sets returned on different processors if this is a distributed vector will form disjoint sets that add up to the complete index set. Obviously, if a vector is created on only one processor, then the result would satisfy

vec.locally_owned_elements() == complete_index_set (vec.size())
IndexSet complete_index_set(const IndexSet::size_type N)
Definition index_set.h:1089

Since the current data type does not support parallel data storage across different processors, the returned index set is the complete index set.

◆ size()

template<typename Number >
size_type Vector< Number >::size ( ) const

Return dimension of the vector.

◆ locally_owned_size()

template<typename Number >
size_type Vector< Number >::locally_owned_size ( ) const

Return local dimension of the vector. Since this vector does not support distributed data this is always the same value as size().

Note
This function exists for compatibility with LinearAlgebra::ReadWriteVector.

◆ all_zero()

template<typename Number >
bool Vector< Number >::all_zero ( ) const

Return whether the vector contains only elements with value zero. This function is mainly for internal consistency checks and should seldom be used when not in debug mode since it uses quite some time.

◆ is_non_negative()

template<typename Number >
bool Vector< Number >::is_non_negative ( ) const

Return true if the vector has no negative entries, i.e. all entries are zero or positive. This function is used, for example, to check whether refinement indicators are really all positive (or zero).

The function obviously only makes sense if the template argument of this class is a real type. If it is a complex type, then an exception is thrown.

◆ memory_consumption()

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

Determine an estimate for the memory consumption (in bytes) of this object.

◆ has_ghost_elements()

template<typename Number >
bool Vector< Number >::has_ghost_elements ( ) const

This function exists for compatibility with the parallel vector classes (e.g., LinearAlgebra::distributed::Vector class) and always returns false since this implementation is serial.

◆ zero_out_ghost_values()

template<typename Number >
void Vector< Number >::zero_out_ghost_values ( ) const

This function exists for compatibility with the parallel vector classes (e.g., LinearAlgebra::distributed::Vector class) and does nothing since this implementation is serial.

◆ maybe_reset_thread_partitioner()

template<typename Number >
void Vector< Number >::maybe_reset_thread_partitioner ( )
private

Convenience function used at the end of initialization or reinitialization. Resets (if necessary) the loop partitioner to the correct state, based on its current state and the length of the vector.

◆ do_reinit()

template<typename Number >
void Vector< Number >::do_reinit ( const size_type  new_size,
const bool  omit_zeroing_entries,
const bool  reset_partitioner 
)
private

Actual implementation of the reinit functions.

◆ lp_norm() [2/2]

Vector< int >::real_type Vector< int >::lp_norm ( const real_type  ) const

Definition at line 80 of file vector.cc.

◆ 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 136 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 156 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 204 of file subscriptor.cc.

◆ 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 53 of file subscriptor.cc.

Friends And Related Symbol Documentation

◆ Vector

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

Typedef for the vector type used

Typedef for the vector type used.

Definition at line 1049 of file vector.h.

◆ swap() [1/3]

template<typename Number >
void swap ( LinearAlgebra::CUDAWrappers::Vector< Number > &  u,
LinearAlgebra::CUDAWrappers::Vector< Number > &  v 
)
related

Global function swap which overloads the default implementation of the C++ standard library which uses a temporary object. The function simply exchanges the data of the two vectors.

Definition at line 373 of file cuda_vector.h.

◆ swap() [2/3]

template<typename Number , typename MemorySpace >
void swap ( LinearAlgebra::distributed::Vector< Number, MemorySpace > &  u,
LinearAlgebra::distributed::Vector< Number, MemorySpace > &  v 
)
related

Global function swap which overloads the default implementation of the C++ standard library which uses a temporary object. The function simply exchanges the data of the two vectors.

Definition at line 1809 of file la_parallel_vector.h.

◆ swap() [3/3]

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

Global function swap which overloads the default implementation of the C++ standard library which uses a temporary object. The function simply exchanges the data of the two vectors.

Definition at line 1225 of file read_write_vector.h.

Member Data Documentation

◆ values

template<typename Number >
AlignedVector<Number> Vector< Number >::values
private

Array of elements owned by this vector.

Definition at line 1022 of file vector.h.

◆ thread_loop_partitioner

template<typename Number >
std::shared_ptr<parallel::internal::TBBPartitioner> Vector< Number >::thread_loop_partitioner
mutableprivate

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

Definition at line 1045 of file 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 219 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 225 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 241 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 249 of file subscriptor.h.

◆ mutex

std::mutex Subscriptor::mutex
staticprivateinherited

A mutex used to ensure data consistency when printing out the list of subscribers.

Definition at line 271 of file subscriptor.h.


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