Reference documentation for deal.II version 9.5.0
|
#include <deal.II/lac/cuda_vector.h>
Public Types | |
using | value_type = typename VectorSpaceVector< Number >::value_type |
using | size_type = typename VectorSpaceVector< Number >::size_type |
using | real_type = typename VectorSpaceVector< Number >::real_type |
Public Member Functions | |
Vector () | |
Vector (const Vector< Number > &V) | |
Vector (Vector< Number > &&) noexcept=default | |
Vector (const size_type n) | |
Vector & | operator= (const Vector< Number > &v) |
Vector & | operator= (Vector< Number > &&v) noexcept=default |
virtual void | swap (Vector< Number > &v) |
void | reinit (const size_type n, const bool omit_zeroing_entries=false) |
virtual void | reinit (const VectorSpaceVector< Number > &V, const bool omit_zeroing_entries=false) override |
virtual void | import_elements (const ReadWriteVector< Number > &V, VectorOperation::values operation, std::shared_ptr< const Utilities::MPI::CommunicationPatternBase > communication_pattern={}) override |
virtual void | import (const ReadWriteVector< Number > &V, VectorOperation::values operation, std::shared_ptr< const Utilities::MPI::CommunicationPatternBase > communication_pattern={}) override |
virtual Vector< Number > & | operator= (const Number s) override |
virtual Vector< Number > & | operator*= (const Number factor) override |
virtual Vector< Number > & | operator/= (const Number factor) override |
virtual Vector< Number > & | operator+= (const VectorSpaceVector< Number > &V) override |
virtual Vector< Number > & | operator-= (const VectorSpaceVector< Number > &V) override |
virtual Number | operator* (const VectorSpaceVector< Number > &V) const override |
virtual void | add (const Number a) override |
virtual void | add (const Number a, const VectorSpaceVector< Number > &V) override |
virtual void | add (const Number a, const VectorSpaceVector< Number > &V, const Number b, const VectorSpaceVector< Number > &W) override |
virtual void | sadd (const Number s, const Number a, const VectorSpaceVector< Number > &V) override |
virtual void | scale (const VectorSpaceVector< Number > &scaling_factors) override |
virtual void | equ (const Number a, const VectorSpaceVector< Number > &V) override |
virtual bool | all_zero () const override |
virtual value_type | mean_value () const override |
virtual real_type | l1_norm () const override |
virtual real_type | l2_norm () const override |
real_type | norm_sqr () const |
virtual real_type | linfty_norm () const override |
virtual Number | add_and_dot (const Number a, const VectorSpaceVector< Number > &V, const VectorSpaceVector< Number > &W) override |
Number * | get_values () const |
virtual size_type | size () const override |
virtual ::IndexSet | locally_owned_elements () const override |
virtual void | print (std::ostream &out, const unsigned int precision=2, const bool scientific=true, const bool across=true) const override |
virtual std::size_t | memory_consumption () const override |
virtual void | compress (VectorOperation::values) |
Static Public Member Functions | |
static ::ExceptionBase & | ExcVectorTypeNotCompatible () |
Private Attributes | |
std::unique_ptr< Number[], void(*)(Number *)> | val |
size_type | n_elements |
This class implements a vector using CUDA for use on Nvidia GPUs. This class is derived from the LinearAlgebra::VectorSpaceVector class.
Definition at line 55 of file cuda_vector.h.
using LinearAlgebra::CUDAWrappers::Vector< Number >::value_type = typename VectorSpaceVector<Number>::value_type |
Definition at line 58 of file cuda_vector.h.
using LinearAlgebra::CUDAWrappers::Vector< Number >::size_type = typename VectorSpaceVector<Number>::size_type |
Definition at line 59 of file cuda_vector.h.
using LinearAlgebra::CUDAWrappers::Vector< Number >::real_type = typename VectorSpaceVector<Number>::real_type |
Definition at line 60 of file cuda_vector.h.
Constructor. Create a vector of dimension zero.
Definition at line 41 of file cuda_vector.cc.
Copy constructor.
Definition at line 49 of file cuda_vector.cc.
|
defaultnoexcept |
Move constructor.
Constructor. Set dimension to n
and initialize all elements with zero.
The constructor is made explicit to avoid accident like this: v=0;
. Presumably, the user wants to set every elements 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.
Definition at line 86 of file cuda_vector.cc.
Vector< Number > & Vector< Number >::operator= | ( | const Vector< Number > & | v | ) |
Copy assignment operator.
Definition at line 66 of file cuda_vector.cc.
|
defaultnoexcept |
Move assignment operator.
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.
Definition at line 411 of file cuda_vector.h.
void Vector< Number >::reinit | ( | const size_type | n, |
const bool | omit_zeroing_entries = false |
||
) |
Reinit functionality. The flag omit_zeroing_entries
determines whether the vector should be filled with zero (false) or left untouched (true).
Definition at line 97 of file cuda_vector.cc.
|
overridevirtual |
Change the dimension to that of the vector V. The elements of V are not copied.
Implements LinearAlgebra::VectorSpaceVector< Number >.
Definition at line 119 of file cuda_vector.cc.
|
overridevirtual |
Import all the element from the input vector V
. VectorOperation::values operation
is used to decide if the elements int V
should be added to the current vector or replace the current elements. The last parameter is not used. It is only used for distributed vectors. This is the function that should be used to copy a vector to the GPU.
Implements LinearAlgebra::VectorSpaceVector< Number >.
Definition at line 129 of file cuda_vector.cc.
|
inlineoverridevirtual |
Implements LinearAlgebra::VectorSpaceVector< Number >.
Definition at line 154 of file cuda_vector.h.
|
overridevirtual |
Sets all elements of the vector to the scalar s
. This operation is only allowed if s
is equal to zero.
Implements LinearAlgebra::VectorSpaceVector< Number >.
Definition at line 175 of file cuda_vector.cc.
|
overridevirtual |
Multiply the entive vector by a fixed factor.
Implements LinearAlgebra::VectorSpaceVector< Number >.
Definition at line 191 of file cuda_vector.cc.
|
overridevirtual |
Divide the entire vector by a fixed factor.
Implements LinearAlgebra::VectorSpaceVector< Number >.
Definition at line 206 of file cuda_vector.cc.
|
overridevirtual |
Add the vector V
to the present one.
Implements LinearAlgebra::VectorSpaceVector< Number >.
Definition at line 222 of file cuda_vector.cc.
|
overridevirtual |
Subtract the vector V
from the present one.
Implements LinearAlgebra::VectorSpaceVector< Number >.
Definition at line 247 of file cuda_vector.cc.
|
overridevirtual |
Return the scalar product of two vectors.
Implements LinearAlgebra::VectorSpaceVector< Number >.
Definition at line 272 of file cuda_vector.cc.
|
overridevirtual |
Add to
all components. Note that a
is a scalar not a vector.
Implements LinearAlgebra::VectorSpaceVector< Number >.
Definition at line 315 of file cuda_vector.cc.
|
overridevirtual |
Simple addition of a multiple of a vector, i.e. *this += a*V
.
Implements LinearAlgebra::VectorSpaceVector< Number >.
Definition at line 328 of file cuda_vector.cc.
|
overridevirtual |
Multiple additions of scaled vectors, i.e. *this += a*V+b*W
.
Implements LinearAlgebra::VectorSpaceVector< Number >.
Definition at line 352 of file cuda_vector.cc.
|
overridevirtual |
Scaling and simple addition of a multiple of a vector, i.e. this = s(*this)+a*V
Implements LinearAlgebra::VectorSpaceVector< Number >.
Definition at line 390 of file cuda_vector.cc.
|
overridevirtual |
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.
Implements LinearAlgebra::VectorSpaceVector< Number >.
Definition at line 417 of file cuda_vector.cc.
|
overridevirtual |
Assignment *this = a*V
.
Implements LinearAlgebra::VectorSpaceVector< Number >.
Definition at line 440 of file cuda_vector.cc.
Return whether the vector contains only elements with value zero.
Implements LinearAlgebra::VectorSpaceVector< Number >.
Definition at line 465 of file cuda_vector.cc.
|
overridevirtual |
Return the mean value of all the entries of this vector.
Implements LinearAlgebra::VectorSpaceVector< Number >.
Definition at line 474 of file cuda_vector.cc.
|
overridevirtual |
Return the l1 norm of the vector (i.e., the sum of the absolute values of all entries among all processors).
Implements LinearAlgebra::VectorSpaceVector< Number >.
Definition at line 505 of file cuda_vector.cc.
|
overridevirtual |
Return the l2 norm of the vector (i.e., the square root of the sum of the square of all entries among all processors).
Implements LinearAlgebra::VectorSpaceVector< Number >.
Definition at line 535 of file cuda_vector.cc.
Vector< Number >::real_type Vector< Number >::norm_sqr |
Return the square of the \(l_2\)-norm.
Definition at line 544 of file cuda_vector.cc.
|
overridevirtual |
Return the maximum norm of the vector (i.e., the maximum absolute value among all entries and among all processors).
Implements LinearAlgebra::VectorSpaceVector< Number >.
Definition at line 553 of file cuda_vector.cc.
|
overridevirtual |
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}\).
Implements LinearAlgebra::VectorSpaceVector< Number >.
Definition at line 583 of file cuda_vector.cc.
|
inline |
Return the pointer to the underlying array. Ownership still resides with this class.
Definition at line 385 of file cuda_vector.h.
|
inlineoverridevirtual |
Return the size of the vector.
Implements LinearAlgebra::VectorSpaceVector< Number >.
Definition at line 394 of file cuda_vector.h.
Return an index set that describe which elements of this vector are owned by the current processor, i.e. [0, size).
Implements LinearAlgebra::VectorSpaceVector< Number >.
Definition at line 402 of file cuda_vector.h.
|
overridevirtual |
Print the vector to the output stream out
.
Implements LinearAlgebra::VectorSpaceVector< Number >.
Definition at line 632 of file cuda_vector.cc.
|
overridevirtual |
Return the memory consumption of this class in bytes.
Implements LinearAlgebra::VectorSpaceVector< Number >.
Definition at line 668 of file cuda_vector.cc.
|
inlinevirtualinherited |
This function does nothing and only exists for backward compatibility.
Reimplemented in LinearAlgebra::distributed::BlockVector< Number >, LinearAlgebra::distributed::Vector< Number, MemorySpace >, and LinearAlgebra::distributed::Vector< double >.
Definition at line 238 of file vector_space_vector.h.
|
private |
Pointer to the array of elements of this vector.
Definition at line 352 of file cuda_vector.h.
|
private |
Number of elements in the vector.
Definition at line 357 of file cuda_vector.h.