Reference documentation for deal.II version 9.0.0
|
#include <deal.II/lac/vector_space_vector.h>
Public Member Functions | |
virtual void | reinit (const VectorSpaceVector< Number > &V, const bool omit_zeroing_entries=false)=0 |
virtual VectorSpaceVector< Number > & | operator= (const Number s)=0 |
virtual VectorSpaceVector< Number > & | operator*= (const Number factor)=0 |
virtual VectorSpaceVector< Number > & | operator/= (const Number factor)=0 |
virtual VectorSpaceVector< Number > & | operator+= (const VectorSpaceVector< Number > &V)=0 |
virtual VectorSpaceVector< Number > & | operator-= (const VectorSpaceVector< Number > &V)=0 |
virtual void | import (const ReadWriteVector< Number > &V, VectorOperation::values operation, std::shared_ptr< const CommunicationPatternBase > communication_pattern=std::shared_ptr< const CommunicationPatternBase >())=0 |
virtual Number | operator* (const VectorSpaceVector< Number > &V) const =0 |
virtual void | add (const Number a)=0 |
virtual void | add (const Number a, const VectorSpaceVector< Number > &V)=0 |
virtual void | add (const Number a, const VectorSpaceVector< Number > &V, const Number b, const VectorSpaceVector< Number > &W)=0 |
virtual void | sadd (const Number s, const Number a, const VectorSpaceVector< Number > &V)=0 |
virtual void | scale (const VectorSpaceVector< Number > &scaling_factors)=0 |
virtual void | equ (const Number a, const VectorSpaceVector< Number > &V)=0 |
virtual bool | all_zero () const =0 |
virtual value_type | mean_value () const =0 |
virtual real_type | l1_norm () const =0 |
virtual real_type | l2_norm () const =0 |
virtual real_type | linfty_norm () const =0 |
virtual Number | add_and_dot (const Number a, const VectorSpaceVector< Number > &V, const VectorSpaceVector< Number > &W)=0 |
virtual void | compress (VectorOperation::values) |
virtual size_type | size () const =0 |
virtual ::IndexSet | locally_owned_elements () const =0 |
virtual void | print (std::ostream &out, const unsigned int precision=3, const bool scientific=true, const bool across=true) const =0 |
virtual std::size_t | memory_consumption () const =0 |
virtual | ~VectorSpaceVector ()=default |
VectorSpaceVector is an abstract class that is used to define the interface that vector classes need to implement when they want to implement global operations. This class is complementary of ReadWriteVector which allows the access of individual elements but does not allow global operations.
Definition at line 51 of file vector_space_vector.h.
|
virtualdefault |
Destructor. Declared as virtual so that inheriting classes (which may manage their own memory) are destroyed correctly.
|
pure virtual |
Change the dimension to that of the vector V. The elements of V are not copied.
Implemented in LinearAlgebra::distributed::Vector< Number >, LinearAlgebra::distributed::Vector< double >, LinearAlgebra::distributed::BlockVector< Number >, LinearAlgebra::Vector< Number >, LinearAlgebra::CUDAWrappers::Vector< Number >, and LinearAlgebra::EpetraWrappers::Vector.
|
pure virtual |
Sets all elements of the vector to the scalar s
. This operation is only allowed if s
is equal to zero.
Implemented in LinearAlgebra::distributed::Vector< Number >, LinearAlgebra::distributed::Vector< double >, LinearAlgebra::Vector< Number >, LinearAlgebra::CUDAWrappers::Vector< Number >, and LinearAlgebra::EpetraWrappers::Vector.
|
pure virtual |
Multiply the entire vector by a fixed factor.
Implemented in LinearAlgebra::distributed::Vector< Number >, LinearAlgebra::distributed::Vector< double >, LinearAlgebra::distributed::BlockVector< Number >, LinearAlgebra::Vector< Number >, LinearAlgebra::EpetraWrappers::Vector, and LinearAlgebra::CUDAWrappers::Vector< Number >.
|
pure virtual |
Divide the entire vector by a fixed factor.
Implemented in LinearAlgebra::distributed::Vector< Number >, LinearAlgebra::distributed::Vector< double >, LinearAlgebra::distributed::BlockVector< Number >, LinearAlgebra::Vector< Number >, LinearAlgebra::EpetraWrappers::Vector, and LinearAlgebra::CUDAWrappers::Vector< Number >.
|
pure virtual |
Add the vector V
to the present one.
Implemented in LinearAlgebra::distributed::Vector< Number >, LinearAlgebra::distributed::Vector< double >, LinearAlgebra::distributed::BlockVector< Number >, LinearAlgebra::Vector< Number >, LinearAlgebra::EpetraWrappers::Vector, and LinearAlgebra::CUDAWrappers::Vector< Number >.
|
pure virtual |
Subtract the vector V
from the present one.
Implemented in LinearAlgebra::distributed::Vector< Number >, LinearAlgebra::distributed::Vector< double >, LinearAlgebra::distributed::BlockVector< Number >, LinearAlgebra::Vector< Number >, LinearAlgebra::EpetraWrappers::Vector, and LinearAlgebra::CUDAWrappers::Vector< Number >.
|
pure virtual |
Import all the elements present in the vector's IndexSet from the input vector V
. 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.
Implemented in LinearAlgebra::distributed::Vector< Number >, LinearAlgebra::distributed::Vector< double >, LinearAlgebra::distributed::BlockVector< Number >, LinearAlgebra::Vector< Number >, LinearAlgebra::EpetraWrappers::Vector, and LinearAlgebra::CUDAWrappers::Vector< Number >.
|
pure virtual |
Return the scalar product of two vectors.
Implemented in LinearAlgebra::distributed::Vector< Number >, LinearAlgebra::distributed::Vector< double >, LinearAlgebra::distributed::BlockVector< Number >, LinearAlgebra::Vector< Number >, LinearAlgebra::EpetraWrappers::Vector, and LinearAlgebra::CUDAWrappers::Vector< Number >.
|
pure virtual |
Add a
to all components. Note that a
is a scalar not a vector.
Implemented in LinearAlgebra::distributed::Vector< Number >, LinearAlgebra::distributed::Vector< double >, LinearAlgebra::distributed::BlockVector< Number >, LinearAlgebra::Vector< Number >, LinearAlgebra::EpetraWrappers::Vector, and LinearAlgebra::CUDAWrappers::Vector< Number >.
|
pure virtual |
Simple addition of a multiple of a vector, i.e. *this += a*V
.
Implemented in LinearAlgebra::distributed::Vector< Number >, LinearAlgebra::distributed::Vector< double >, LinearAlgebra::distributed::BlockVector< Number >, LinearAlgebra::Vector< Number >, LinearAlgebra::EpetraWrappers::Vector, and LinearAlgebra::CUDAWrappers::Vector< Number >.
|
pure virtual |
Multiple addition of scaled vectors, i.e. *this += a*V+b*W
.
Implemented in LinearAlgebra::distributed::Vector< Number >, LinearAlgebra::distributed::Vector< double >, LinearAlgebra::distributed::BlockVector< Number >, LinearAlgebra::Vector< Number >, LinearAlgebra::EpetraWrappers::Vector, and LinearAlgebra::CUDAWrappers::Vector< Number >.
|
pure virtual |
Scaling and simple addition of a multiple of a vector, i.e. *this = s*(*this)+a*V
.
Implemented in LinearAlgebra::distributed::Vector< Number >, LinearAlgebra::distributed::Vector< double >, LinearAlgebra::distributed::BlockVector< Number >, LinearAlgebra::Vector< Number >, LinearAlgebra::EpetraWrappers::Vector, and LinearAlgebra::CUDAWrappers::Vector< Number >.
|
pure virtual |
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.
Implemented in LinearAlgebra::distributed::Vector< Number >, LinearAlgebra::distributed::Vector< double >, LinearAlgebra::distributed::BlockVector< Number >, LinearAlgebra::Vector< Number >, LinearAlgebra::EpetraWrappers::Vector, and LinearAlgebra::CUDAWrappers::Vector< Number >.
|
pure virtual |
Assignment *this = a*V
.
Implemented in LinearAlgebra::distributed::Vector< Number >, LinearAlgebra::distributed::Vector< double >, LinearAlgebra::distributed::BlockVector< Number >, LinearAlgebra::Vector< Number >, LinearAlgebra::EpetraWrappers::Vector, and LinearAlgebra::CUDAWrappers::Vector< Number >.
|
pure virtual |
Return whether the vector contains only elements with value zero.
Implemented in LinearAlgebra::distributed::Vector< Number >, LinearAlgebra::distributed::Vector< double >, LinearAlgebra::distributed::BlockVector< Number >, LinearAlgebra::Vector< Number >, LinearAlgebra::EpetraWrappers::Vector, and LinearAlgebra::CUDAWrappers::Vector< Number >.
|
pure virtual |
Return the mean value of all the entries of this vector.
Implemented in LinearAlgebra::distributed::Vector< Number >, LinearAlgebra::distributed::Vector< double >, LinearAlgebra::distributed::BlockVector< Number >, LinearAlgebra::Vector< Number >, LinearAlgebra::EpetraWrappers::Vector, and LinearAlgebra::CUDAWrappers::Vector< Number >.
|
pure virtual |
Return the l1 norm of the vector (i.e., the sum of the absolute values of all entries among all processors).
Implemented in LinearAlgebra::distributed::Vector< Number >, LinearAlgebra::distributed::Vector< double >, LinearAlgebra::distributed::BlockVector< Number >, LinearAlgebra::Vector< Number >, LinearAlgebra::EpetraWrappers::Vector, and LinearAlgebra::CUDAWrappers::Vector< Number >.
|
pure virtual |
Return the l2 norm of the vector (i.e., the square root of the sum of the square of all entries among all processors).
Implemented in LinearAlgebra::distributed::Vector< Number >, LinearAlgebra::distributed::Vector< double >, LinearAlgebra::distributed::BlockVector< Number >, LinearAlgebra::Vector< Number >, LinearAlgebra::EpetraWrappers::Vector, and LinearAlgebra::CUDAWrappers::Vector< Number >.
|
pure virtual |
Return the maximum norm of the vector (i.e., the maximum absolute value among all entries and among all processors).
Implemented in LinearAlgebra::distributed::Vector< Number >, LinearAlgebra::distributed::Vector< double >, LinearAlgebra::distributed::BlockVector< Number >, LinearAlgebra::Vector< Number >, LinearAlgebra::EpetraWrappers::Vector, and LinearAlgebra::CUDAWrappers::Vector< Number >.
|
pure virtual |
Perform 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}\).
Implemented in LinearAlgebra::distributed::Vector< Number >, LinearAlgebra::distributed::Vector< double >, LinearAlgebra::distributed::BlockVector< Number >, LinearAlgebra::Vector< Number >, LinearAlgebra::EpetraWrappers::Vector, and LinearAlgebra::CUDAWrappers::Vector< Number >.
|
inlinevirtual |
This function does nothing and only exists for backward compatibility.
Definition at line 198 of file vector_space_vector.h.
|
pure virtual |
Return the global size of the vector, equal to the sum of the number of locally owned indices among all processors.
Implemented in LinearAlgebra::distributed::Vector< Number >, LinearAlgebra::distributed::Vector< double >, LinearAlgebra::distributed::BlockVector< Number >, LinearAlgebra::Vector< Number >, LinearAlgebra::EpetraWrappers::Vector, and LinearAlgebra::CUDAWrappers::Vector< Number >.
|
pure virtual |
Return an index set that describes which elements of this vector are owned by the current 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
Implemented in LinearAlgebra::distributed::Vector< Number >, LinearAlgebra::distributed::Vector< double >, LinearAlgebra::distributed::BlockVector< Number >, LinearAlgebra::Vector< Number >, LinearAlgebra::EpetraWrappers::Vector, and LinearAlgebra::CUDAWrappers::Vector< Number >.
|
pure virtual |
Print the vector to the output stream out
.
Implemented in LinearAlgebra::distributed::Vector< Number >, LinearAlgebra::distributed::Vector< double >, LinearAlgebra::distributed::BlockVector< Number >, LinearAlgebra::Vector< Number >, LinearAlgebra::EpetraWrappers::Vector, and LinearAlgebra::CUDAWrappers::Vector< Number >.
|
pure virtual |
Return the memory consumption of this class in bytes.
Implemented in LinearAlgebra::distributed::Vector< Number >, LinearAlgebra::distributed::Vector< double >, LinearAlgebra::distributed::BlockVector< Number >, LinearAlgebra::Vector< Number >, LinearAlgebra::EpetraWrappers::Vector, and LinearAlgebra::CUDAWrappers::Vector< Number >.