Reference documentation for deal.II version 9.1.1
|
#include <deal.II/lac/vector.h>
Public Types | |
using | value_type = Number |
using | real_type = typename numbers::NumberTraits< Number >::real_type |
Public Member Functions | |
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 | ratio (const Vector< Number > &a, const Vector< Number > &b) |
void | update_ghost_values () const |
Input and output | |
void | print (const char *format=nullptr) const |
void | print (std::ostream &out, const unsigned int precision=3, const bool scientific=true, const bool across=true) const |
void | print (LogStream &out, const unsigned int width=6, 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) |
Information about the object | |
bool | in_local_range (const size_type global_index) const |
IndexSet | locally_owned_elements () const |
size_type | size () const |
bool | all_zero () const |
bool | is_non_negative () const |
std::size_t | memory_consumption () const |
Public Member Functions inherited from Subscriptor | |
Subscriptor () | |
Subscriptor (const Subscriptor &) | |
Subscriptor (Subscriptor &&) noexcept | |
virtual | ~Subscriptor () |
Subscriptor & | operator= (const Subscriptor &) |
Subscriptor & | operator= (Subscriptor &&) noexcept |
void | subscribe (std::atomic< bool > *const validity, const std::string &identifier="") const |
void | unsubscribe (std::atomic< bool > *const validity, const std::string &identifier="") const |
unsigned int | n_subscriptions () const |
template<typename StreamType > | |
void | list_subscribers (StreamType &stream) const |
void | list_subscribers () const |
template<class Archive > | |
void | serialize (Archive &ar, const unsigned int version) |
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) |
Private Attributes | |
AlignedVector< Number > | values |
std::shared_ptr< parallel::internal::TBBPartitioner > | thread_loop_partitioner |
Friends | |
template<typename Number2 > | |
class | Vector |
Related Functions | |
(Note that these are not member functions.) | |
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) |
Additional Inherited Members | |
Static Public Member Functions inherited from Subscriptor | |
static ::ExceptionBase & | ExcInUse (int arg1, std::string arg2, std::string arg3) |
static ::ExceptionBase & | ExcNoSubscriber (std::string arg1, std::string arg2) |
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.
<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 36 of file function.h.
using Vector< Number >::value_type = 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.
Constructor. Create a vector of dimension zero.
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.
Move constructor. Creates a new vector by stealing the internal data of the vector 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.
|
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:
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.
|
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.
|
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.
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
. 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< 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.
Destructor, deallocates memory. Made virtual to allow for derived classes to behave properly.
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.
|
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 behaviour 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.
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 \]
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.
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 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.
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.
Copy the given vector. Resize the present vector if necessary.
|
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.
Vector<Number>& Vector< Number >::operator= | ( | const Vector< Number2 > & | v | ) |
Copy the given vector. Resize the present vector if necessary.
Vector<Number>& Vector< Number >::operator= | ( | const BlockVector< Number > & | v | ) |
Copy operator for assigning a block vector to a regular vector.
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.
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.
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
. 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.
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.
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}\).
Return the square of the \(l_2\)-norm.
Number Vector< Number >::mean_value | ( | ) | const |
Mean value of the elements of this vector.
\(l_1\)-norm of the vector. The sum of the absolute values.
\(l_2\)-norm of the vector. The square root of the sum of the squares of the elements.
\(l_p\)-norm of the vector. The pth root of the sum of the pth powers of the absolute values of the elements.
Maximum absolute value of the elements.
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
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}\).
pointer Vector< Number >::data | ( | ) |
Return a pointer to the underlying data buffer.
const_pointer Vector< Number >::data | ( | ) | const |
Return a const pointer to the underlying data buffer.
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.
const_iterator Vector< Number >::begin | ( | ) | const |
Return constant iterator to the start of the vectors.
iterator Vector< Number >::end | ( | ) |
Return an iterator pointing to the element past the end of the array.
const_iterator Vector< Number >::end | ( | ) | const |
Return a constant iterator pointing to the element past the end of the array.
Number Vector< Number >::operator() | ( | const size_type | i | ) | const |
Access the value of the ith
component.
Number& Vector< Number >::operator() | ( | const size_type | i | ) |
Access the ith
component as a writeable reference.
Number Vector< Number >::operator[] | ( | const size_type | i | ) | const |
Access the value of the ith
component.
Exactly the same as operator().
Number& Vector< Number >::operator[] | ( | const size_type | i | ) |
Access the ith
component as a writeable reference.
Exactly the same as operator().
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
indices
and values
arrays must be identical. 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
values_begin
as there are iterators between indices_begin
and indices_end
. Vector<Number>& Vector< Number >::operator+= | ( | const Vector< Number > & | V | ) |
Add the given vector to the present one.
Vector<Number>& Vector< Number >::operator-= | ( | const Vector< Number > & | V | ) |
Subtract the given vector from the present one.
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
.
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.
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.
void Vector< Number >::add | ( | const Number | s | ) |
Addition of s
to all components. Note that s
is a scalar and not a vector.
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
.
Simple addition of a multiple of a vector, i.e. *this += a*V
.
void Vector< Number >::sadd | ( | const Number | s, |
const Vector< Number > & | V | ||
) |
Scaling and simple vector addition, i.e. *this = s*(*this)+V
.
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
.
Scale each element of the vector by a constant value.
Scale each element of the vector by the inverse of the given value.
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.
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.
Assignment *this = a*u
.
void Vector< Number >::equ | ( | const Number | a, |
const Vector< Number2 > & | u | ||
) |
Assignment *this = a*u
.
void Vector< Number >::ratio | ( | const Vector< Number > & | a, |
const Vector< Number > & | b | ||
) |
Compute the elementwise ratio of the two given vectors, that is let this[i] = a[i]/b[i]
. This is useful for example if you want to compute the cellwise ratio of true to estimated error.
This vector is appropriately scaled to hold the result.
If any of the b[i]
is zero, the result is undefined. No attempt is made to catch such situations.
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).
void Vector< Number >::print | ( | const char * | format = nullptr | ) | const |
Output of vector in user-defined format. For complex-valued vectors, the format should include specifiers for both the real and imaginary parts.
This function is deprecated.
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.
void Vector< Number >::print | ( | LogStream & | out, |
const unsigned int | width = 6 , |
||
const bool | across = true |
||
) | const |
Print to a LogStream. width
is used as argument to the std::setw manipulator, if printing across. 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.
This function is deprecated.
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.
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.
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.
void Vector< Number >::load | ( | Archive & | ar, |
const unsigned int | version | ||
) |
Read the data of this object from a stream for the purpose of serialization.
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.
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
Since the current data type does not support parallel data storage across different processors, the returned index set is the complete index set.
size_type Vector< Number >::size | ( | ) | const |
Return dimension of the vector.
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.
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.
std::size_t Vector< Number >::memory_consumption | ( | ) | const |
Determine an estimate for the memory consumption (in bytes) of this object.
|
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.
|
private |
Actual implementation of the reinit functions.
Make all other vector types friends.
Typedef for the vector type used
Typedef for the vector type used.
Implementation of a parallel vector class. The design of this class is similar to the standard Vector class in deal.II, with the exception that storage is distributed with MPI.
The vector is designed for the following scheme of parallel partitioning:
[my_first_index,my_last_index)
. reinit (locally_owned, ghost_indices, communicator)
, and retained until the partitioning is changed. This allows for efficient parallel communication of indices. In particular, it stores the communication pattern, rather than having to compute it again for every communication. For more information on ghost vectors, see also the glossary entry on vectors with ghost elements. local_element()
. Locally owned indices are placed first, [0, local_size()), and then all ghost indices follow after them contiguously, [local_size(), local_size()+n_ghost_entries()). Functions related to parallel functionality:
compress()
goes through the data associated with ghost indices and communicates it to the owner process, which can then add it to the correct position. This can be used e.g. after having run an assembly routine involving ghosts that fill this vector. Note that the insert
mode of compress()
does not set the elements included in ghost entries but simply discards them, assuming that the owning processor has set them to the desired value already (See also the glossary entry on compress). update_ghost_values()
function imports the data from the owning processor to the ghost indices in order to provide read access to the data associated with ghosts. This vector can take two different states with respect to ghost elements:
operator= (0.)
), the vector does only allow writing into ghost elements but not reading from ghost elements. true
exactly when ghost elements have been updated and false
otherwise, irrespective of the actual number of ghost entries in the vector layout (for that information, use n_ghost_entries() instead). This vector uses the facilities of the class Vector<Number> for implementing the operations on the local range of the vector. In particular, it also inherits thread parallelism that splits most vector-vector operations into smaller chunks if the program uses multiple threads. This may or may not be desired when working also with MPI.
This vector class is based on two different number types for indexing. The so-called global index type encodes the overall size of the vector. Its type is types::global_dof_index. The largest possible value is 2^32-1
or approximately 4 billion in case 64 bit integers are disabled at configuration of deal.II (default case) or 2^64-1
or approximately 10^19
if 64 bit integers are enabled (see the glossary entry on When to use types::global_dof_index instead of unsigned int for further information).
The second relevant index type is the local index used within one MPI rank. As opposed to the global index, the implementation assumes 32-bit unsigned integers unconditionally. In other words, to actually use a vector with more than four billion entries, you need to use MPI with more than one rank (which in general is a safe assumption since four billion entries consume at least 16 GB of memory for floats or 32 GB of memory for doubles) and enable 64-bit indices. If more than 4 billion local elements are present, the implementation tries to detect that, which triggers an exception and aborts the code. Note, however, that the detection of overflow is tricky and the detection mechanism might fail in some circumstances. Therefore, it is strongly recommended to not rely on this class to automatically detect the unsupported case.
|
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 360 of file cuda_vector.h.
|
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 1845 of file la_parallel_vector.h.
|
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 1046 of file read_write_vector.h.
|
private |
|
mutableprivate |