Reference documentation for deal.II version 9.0.0
|
#include <deal.II/lac/la_vector.h>
Public Member Functions | |
Vector ()=default | |
Vector (const Vector< Number > &V) | |
Vector (const size_type n) | |
template<typename InputIterator > | |
Vector (const InputIterator first, const InputIterator last) | |
virtual void | reinit (const size_type size, const bool omit_zeroing_entries=false) override |
template<typename Number2 > | |
void | reinit (const ReadWriteVector< Number2 > &in_vector, const bool omit_zeroing_entries=false) |
virtual void | reinit (const IndexSet &locally_stored_indices, const bool omit_zeroing_entries=false) override |
virtual void | reinit (const VectorSpaceVector< Number > &V, const bool omit_zeroing_entries=false) override |
Vector< Number > & | operator= (const Vector< Number > &in_vector) |
template<typename Number2 > | |
Vector< Number > & | operator= (const Vector< Number2 > &in_vector) |
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 | import (const ReadWriteVector< Number > &V, VectorOperation::values operation, std::shared_ptr< const CommunicationPatternBase > communication_pattern=std::shared_ptr< const CommunicationPatternBase >()) 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 VectorSpaceVector< Number >::real_type | l1_norm () const override |
virtual VectorSpaceVector< Number >::real_type | l2_norm () const override |
virtual VectorSpaceVector< Number >::real_type | linfty_norm () const override |
virtual Number | add_and_dot (const Number a, const VectorSpaceVector< Number > &V, const VectorSpaceVector< Number > &W) override |
virtual size_type | size () const override |
virtual ::IndexSet | locally_owned_elements () const override |
virtual void | print (std::ostream &out, const unsigned int precision=3, const bool scientific=true, const bool across=true) const override |
void | block_write (std::ostream &out) const |
void | block_read (std::istream &in) |
virtual std::size_t | memory_consumption () const override |
Public Member Functions inherited from LinearAlgebra::ReadWriteVector< Number > | |
ReadWriteVector () | |
ReadWriteVector (const ReadWriteVector< Number > &in_vector) | |
ReadWriteVector (const size_type size) | |
ReadWriteVector (const IndexSet &locally_stored_indices) | |
~ReadWriteVector ()=default | |
template<typename Number2 > | |
void | reinit (const ReadWriteVector< Number2 > &in_vector, const bool omit_zeroing_entries=false) |
void | reinit (const TrilinosWrappers::MPI::Vector &trilinos_vec) |
template<typename Functor > | |
void | apply (const Functor &func) |
void | swap (ReadWriteVector< Number > &v) |
ReadWriteVector< Number > & | operator= (const ReadWriteVector< Number > &in_vector) |
template<typename Number2 > | |
ReadWriteVector< Number > & | operator= (const ReadWriteVector< Number2 > &in_vector) |
ReadWriteVector< Number > & | operator= (const Number s) |
void | import (const distributed::Vector< Number > &vec, VectorOperation::values operation, const std::shared_ptr< const CommunicationPatternBase > &communication_pattern=std::shared_ptr< const CommunicationPatternBase >()) |
void | import (const PETScWrappers::MPI::Vector &petsc_vec, VectorOperation::values operation, const std::shared_ptr< const CommunicationPatternBase > &communication_pattern=std::shared_ptr< const CommunicationPatternBase >()) |
void | import (const TrilinosWrappers::MPI::Vector &trilinos_vec, VectorOperation::values operation, std::shared_ptr< const CommunicationPatternBase > communication_pattern=std::shared_ptr< const CommunicationPatternBase >()) |
void | import (const EpetraWrappers::Vector &epetra_vec, VectorOperation::values operation, std::shared_ptr< const CommunicationPatternBase > communication_pattern=std::shared_ptr< const CommunicationPatternBase >()) |
void | import (const CUDAWrappers::Vector< Number > &cuda_vec, VectorOperation::values operation, std::shared_ptr< const CommunicationPatternBase > communication_pattern=std::shared_ptr< const CommunicationPatternBase >()) |
size_type | size () const |
size_type | n_elements () const |
const IndexSet & | get_stored_elements () const |
iterator | begin () |
const_iterator | begin () const |
iterator | end () |
const_iterator | end () const |
Number | operator() (const size_type global_index) const |
Number & | operator() (const size_type global_index) |
Number | operator[] (const size_type global_index) const |
Number & | operator[] (const size_type global_index) |
template<typename Number2 > | |
void | extract_subvector_to (const std::vector< size_type > &indices, std::vector< Number2 > &values) const |
template<typename ForwardIterator , typename OutputIterator > | |
void | extract_subvector_to (ForwardIterator indices_begin, const ForwardIterator indices_end, OutputIterator values_begin) const |
Number | local_element (const size_type local_index) const |
Number & | local_element (const size_type local_index) |
template<typename Number2 > | |
void | add (const std::vector< size_type > &indices, const std::vector< Number2 > &values) |
template<typename Number2 > | |
void | add (const std::vector< size_type > &indices, const ReadWriteVector< Number2 > &values) |
template<typename Number2 > | |
void | add (const size_type n_elements, const size_type *indices, const Number2 *values) |
void | print (std::ostream &out, const unsigned int precision=3, const bool scientific=true) 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 (const char *identifier=nullptr) const |
void | unsubscribe (const char *identifier=nullptr) const |
unsigned int | n_subscriptions () const |
void | list_subscribers () const |
template<class Archive > | |
void | serialize (Archive &ar, const unsigned int version) |
Public Member Functions inherited from LinearAlgebra::VectorSpaceVector< Number > | |
virtual void | compress (VectorOperation::values) |
virtual | ~VectorSpaceVector ()=default |
Static Public Member Functions | |
static ::ExceptionBase & | ExcVectorTypeNotCompatible () |
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) |
Private Member Functions | |
template<typename Archive > | |
void | serialize (Archive &ar, const unsigned int version) |
Friends | |
template<typename Number2 > | |
class | Vector |
Additional Inherited Members | |
Public Types inherited from LinearAlgebra::ReadWriteVector< Number > | |
typedef Number | value_type |
Protected Member Functions inherited from LinearAlgebra::ReadWriteVector< Number > | |
void | import (const Epetra_MultiVector &multivector, const IndexSet &locally_owned_elements, VectorOperation::values operation, const MPI_Comm &mpi_comm, const std::shared_ptr< const CommunicationPatternBase > &communication_pattern) |
unsigned int | global_to_local (const types::global_dof_index global_index) const |
void | resize_val (const size_type new_allocated_size) |
EpetraWrappers::CommunicationPattern | create_epetra_comm_pattern (const IndexSet &source_index_set, const MPI_Comm &mpi_comm) |
Protected Attributes inherited from LinearAlgebra::ReadWriteVector< Number > | |
IndexSet | stored_elements |
IndexSet | source_stored_elements |
std::shared_ptr< CommunicationPatternBase > | comm_pattern |
std::unique_ptr< Number[], decltype(free) * > | values |
std::shared_ptr< ::parallel::internal::TBBPartitioner > | thread_loop_partitioner |
Numerical vector of data. This class derives from both LinearAlgebra::ReadWriteVector and LinearAlgebra::VectorSpaceVector. As opposed to the array of the C++ standard library, this class implements an element of a vector space suitable for numerical computations.
Definition at line 63 of file la_vector.h.
|
default |
Constructor. Create a vector of dimension zero.
Copy constructor. Sets the dimension to that of the given vector and copies all elements.
Definition at line 361 of file la_vector.h.
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 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.
Definition at line 370 of file la_vector.h.
|
inline |
Initialize the vector with a given range of values pointed to by the iterators. This function exists in analogy to the std::vector
class.
Definition at line 380 of file la_vector.h.
|
overridevirtual |
Set the global size of the vector to size
. The stored elements have their index in [0,size).
If the flag omit_zeroing_entries
is set to false, the memory will be initialized with zero, otherwise the memory will be untouched (and the user must make sure to fill it with reasonable data before using it).
Reimplemented from LinearAlgebra::ReadWriteVector< Number >.
void LinearAlgebra::Vector< Number >::reinit | ( | const ReadWriteVector< Number2 > & | in_vector, |
const bool | omit_zeroing_entries = false |
||
) |
Uses the same IndexSet as the one of the input vector in_vector
and allocates memory for this vector.
If the flag omit_zeroing_entries
is set to false, the memory will be initialized with zero, otherwise the memory will be untouched (and the user must make sure to fill it with reasonable data before using it).
|
overridevirtual |
Initializes the vector. The indices are specified by locally_stored_indices
.
If the flag omit_zeroing_entries
is set to false, the memory will be initialized with zero, otherwise the memory will be untouched (and the user must make sure to fill it with reasonable data before using it). locally_stored_indices.
Reimplemented from LinearAlgebra::ReadWriteVector< Number >.
|
overridevirtual |
Change the dimension to that of the vector V. The elements of V are not copied.
Implements LinearAlgebra::VectorSpaceVector< Number >.
Vector<Number>& LinearAlgebra::Vector< Number >::operator= | ( | const Vector< Number > & | in_vector | ) |
Copies the data of the input vector in_vector
.
Vector<Number>& LinearAlgebra::Vector< Number >::operator= | ( | const Vector< Number2 > & | in_vector | ) |
Copies the data of the input vector in_vector
.
|
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 >.
|
overridevirtual |
Multiply the entire vector by a fixed factor.
Implements LinearAlgebra::VectorSpaceVector< Number >.
|
overridevirtual |
Divide the entire vector by a fixed factor.
Implements LinearAlgebra::VectorSpaceVector< Number >.
|
overridevirtual |
Add the vector V
to the present one.
Implements LinearAlgebra::VectorSpaceVector< Number >.
|
overridevirtual |
Subtract the vector V
from the present one.
Implements LinearAlgebra::VectorSpaceVector< Number >.
|
overridevirtual |
Return the scalar product of two vectors.
Implements LinearAlgebra::VectorSpaceVector< Number >.
|
overridevirtual |
This function is not implemented and will throw an exception.
Implements LinearAlgebra::VectorSpaceVector< Number >.
|
overridevirtual |
Add a
to all components. Note that a
is a scalar not a vector.
Implements LinearAlgebra::VectorSpaceVector< Number >.
|
overridevirtual |
Simple addition of a multiple of a vector, i.e. *this += a*V
.
Implements LinearAlgebra::VectorSpaceVector< Number >.
|
overridevirtual |
Multiple addition of a multiple of a vector, i.e. *this += a*V+b*W
.
Implements LinearAlgebra::VectorSpaceVector< Number >.
|
overridevirtual |
Scaling and simple addition of a multiple of a vector, i.e. *this = s*(*this)+a*V
.
Implements LinearAlgebra::VectorSpaceVector< Number >.
|
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 >.
|
overridevirtual |
Assignment *this = a*V
.
Implements LinearAlgebra::VectorSpaceVector< Number >.
|
overridevirtual |
Return whether the vector contains only elements with value zero.
Implements LinearAlgebra::VectorSpaceVector< Number >.
|
overridevirtual |
Return the mean value of all the entries of this vector.
Implements LinearAlgebra::VectorSpaceVector< Number >.
|
overridevirtual |
Return the l1 norm of the vector (i.e., the sum of the absolute values of all entries).
Implements LinearAlgebra::VectorSpaceVector< Number >.
|
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 >.
|
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 >.
|
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 >.
|
inlineoverridevirtual |
Return the global size of the vector, equal to the sum of the number of locally owned indices among all processors.
Implements LinearAlgebra::VectorSpaceVector< Number >.
Definition at line 390 of file la_vector.h.
|
inlineoverridevirtual |
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
Implements LinearAlgebra::VectorSpaceVector< Number >.
Definition at line 399 of file la_vector.h.
|
inlineoverridevirtual |
Prints the vector to the output stream out
.
Implements LinearAlgebra::VectorSpaceVector< Number >.
Definition at line 408 of file la_vector.h.
void LinearAlgebra::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 LinearAlgebra::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.
|
inlineoverridevirtual |
Return the memory consumption of this class in bytes.
Implements LinearAlgebra::VectorSpaceVector< Number >.
Definition at line 436 of file la_vector.h.
|
inlineprivate |
Serialize the data of this object using boost. This function is necessary to use boost::archive::text_iarchive and boost::archive::text_oarchive.
Definition at line 421 of file la_vector.h.
Make all other ReadWriteVector types friends.
Typedef for the vector type used.
Definition at line 353 of file la_vector.h.