Reference documentation for deal.II version 9.6.0
|
#include <deal.II/lac/cuda_vector.h>
Public Types | |
using | value_type = Number |
using | size_type = types::global_dof_index |
using | real_type = typename numbers::NumberTraits<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 |
void | swap (Vector< Number > &v) |
void | reinit (const size_type n, const bool omit_zeroing_entries=false) |
void | reinit (const Vector< Number > &V, const bool omit_zeroing_entries=false) |
void | import_elements (const ReadWriteVector< Number > &V, const VectorOperation::values operation, const std::shared_ptr< const Utilities::MPI::CommunicationPatternBase > &communication_pattern={}) |
void | import (const ReadWriteVector< Number > &V, VectorOperation::values operation, std::shared_ptr< const Utilities::MPI::CommunicationPatternBase > communication_pattern={}) |
Vector< Number > & | operator= (const Number s) |
Vector< Number > & | operator*= (const Number factor) |
Vector< Number > & | operator/= (const Number factor) |
Vector< Number > & | operator+= (const Vector< Number > &V) |
Vector< Number > & | operator-= (const Vector< Number > &V) |
Number | operator* (const Vector< Number > &V) const |
void | add (const Number a) |
void | add (const Number a, const Vector< Number > &V) |
void | add (const Number a, const Vector< Number > &V, const Number b, const Vector< Number > &W) |
void | sadd (const Number s, const Number a, const Vector< Number > &V) |
void | scale (const Vector< Number > &scaling_factors) |
void | equ (const Number a, const Vector< Number > &V) |
bool | all_zero () const |
value_type | mean_value () const |
real_type | l1_norm () const |
real_type | l2_norm () const |
real_type | norm_sqr () const |
real_type | linfty_norm () const |
Number | add_and_dot (const Number a, const Vector< Number > &V, const Vector< Number > &W) |
Number * | get_values () const |
size_type | size () const |
::IndexSet | locally_owned_elements () const |
void | print (std::ostream &out, const unsigned int precision=2, const bool scientific=true, const bool across=true) const |
std::size_t | memory_consumption () const |
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.
Definition at line 55 of file cuda_vector.h.
using LinearAlgebra::CUDAWrappers::Vector< Number >::value_type = Number |
Definition at line 58 of file cuda_vector.h.
using LinearAlgebra::CUDAWrappers::Vector< Number >::size_type = types::global_dof_index |
Definition at line 59 of file cuda_vector.h.
using LinearAlgebra::CUDAWrappers::Vector< Number >::real_type = typename numbers::NumberTraits<Number>::real_type |
Definition at line 60 of file cuda_vector.h.
Constructor. Create a vector of dimension zero.
Definition at line 40 of file cuda_vector.cc.
Copy constructor.
Definition at line 48 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 85 of file cuda_vector.cc.
Vector< Number > & Vector< Number >::operator= | ( | const Vector< Number > & | v | ) |
Copy assignment operator.
Definition at line 65 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.
Definition at line 405 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 zeros (false
) or left in an undetermined state (true
).
Definition at line 96 of file cuda_vector.cc.
void Vector< Number >::reinit | ( | const Vector< Number > & | V, |
const bool | omit_zeroing_entries = false ) |
Change the dimension to that of the vector V. The elements of V are not copied.
Definition at line 118 of file cuda_vector.cc.
void Vector< Number >::import_elements | ( | const ReadWriteVector< Number > & | V, |
const VectorOperation::values | operation, | ||
const std::shared_ptr< const Utilities::MPI::CommunicationPatternBase > & | communication_pattern = {} ) |
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.
Definition at line 128 of file cuda_vector.cc.
|
inline |
Definition at line 150 of file cuda_vector.h.
Sets all elements of the vector to the scalar s
. This operation is only allowed if s
is equal to zero.
Definition at line 174 of file cuda_vector.cc.
Multiply the entive vector by a fixed factor.
Definition at line 190 of file cuda_vector.cc.
Divide the entire vector by a fixed factor.
Definition at line 205 of file cuda_vector.cc.
Vector< Number > & Vector< Number >::operator+= | ( | const Vector< Number > & | V | ) |
Add the vector V
to the present one.
Definition at line 221 of file cuda_vector.cc.
Vector< Number > & Vector< Number >::operator-= | ( | const Vector< Number > & | V | ) |
Subtract the vector V
from the present one.
Definition at line 240 of file cuda_vector.cc.
Return the scalar product of two vectors.
Definition at line 259 of file cuda_vector.cc.
void Vector< Number >::add | ( | const Number | a | ) |
Add to
all components. Note that a
is a scalar not a vector.
Definition at line 296 of file cuda_vector.cc.
void Vector< Number >::add | ( | const Number | a, |
const Vector< Number > & | V ) |
Simple addition of a multiple of a vector, i.e. *this += a*V
.
Definition at line 309 of file cuda_vector.cc.
void Vector< Number >::add | ( | const Number | a, |
const Vector< Number > & | V, | ||
const Number | b, | ||
const Vector< Number > & | W ) |
Multiple additions of scaled vectors, i.e. *this += a*V+b*W
.
Definition at line 327 of file cuda_vector.cc.
void Vector< Number >::sadd | ( | const Number | s, |
const Number | a, | ||
const Vector< Number > & | V ) |
Scaling and simple addition of a multiple of a vector, i.e. this = s(*this)+a*V
Definition at line 353 of file cuda_vector.cc.
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.
Definition at line 374 of file cuda_vector.cc.
void Vector< Number >::equ | ( | const Number | a, |
const Vector< Number > & | V ) |
Assignment *this = a*V
.
Definition at line 392 of file cuda_vector.cc.
Return whether the vector contains only elements with value zero.
Definition at line 413 of file cuda_vector.cc.
Vector< Number >::value_type Vector< Number >::mean_value | ( | ) | const |
Return the mean value of all the entries of this vector.
Definition at line 422 of file cuda_vector.cc.
Vector< Number >::real_type Vector< Number >::l1_norm | ( | ) | const |
Return the l1 norm of the vector (i.e., the sum of the absolute values of all entries among all processors).
Definition at line 453 of file cuda_vector.cc.
Vector< Number >::real_type Vector< Number >::l2_norm | ( | ) | const |
Return the l2 norm of the vector (i.e., the square root of the sum of the square of all entries among all processors).
Definition at line 483 of file cuda_vector.cc.
Vector< Number >::real_type Vector< Number >::norm_sqr | ( | ) | const |
Return the square of the \(l_2\)-norm.
Definition at line 492 of file cuda_vector.cc.
Vector< Number >::real_type Vector< Number >::linfty_norm | ( | ) | const |
Return the maximum norm of the vector (i.e., the maximum absolute value among all entries and among all processors).
Definition at line 501 of file cuda_vector.cc.
Number Vector< Number >::add_and_dot | ( | const Number | a, |
const Vector< Number > & | V, | ||
const Vector< Number > & | W ) |
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}\).
Definition at line 531 of file cuda_vector.cc.
|
inline |
Return the pointer to the underlying array. Ownership still resides with this class.
Definition at line 379 of file cuda_vector.h.
|
inline |
Return the size of the vector.
Definition at line 388 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).
Definition at line 396 of file cuda_vector.h.
void Vector< Number >::print | ( | std::ostream & | out, |
const unsigned int | precision = 2, | ||
const bool | scientific = true, | ||
const bool | across = true ) const |
Print the vector to the output stream out
.
Definition at line 566 of file cuda_vector.cc.
std::size_t Vector< Number >::memory_consumption | ( | ) | const |
Return the memory consumption of this class in bytes.
Definition at line 602 of file cuda_vector.cc.
|
private |
Pointer to the array of elements of this vector.
Definition at line 346 of file cuda_vector.h.
|
private |
Number of elements in the vector.
Definition at line 351 of file cuda_vector.h.