Reference documentation for deal.II version 8.5.1
Public Member Functions | Static Public Member Functions | Private Member Functions | Friends | List of all members
LinearAlgebra::Vector< Number > Class Template Reference

#include <deal.II/lac/la_vector.h>

Inheritance diagram for LinearAlgebra::Vector< Number >:
[legend]

Public Member Functions

 Vector ()
 
 Vector (const Vector< Number > &V)
 
 Vector (const size_type n)
 
template<typename InputIterator >
 Vector (const InputIterator first, const InputIterator last)
 
Vector< Number > & operator= (const Vector< Number > &in_vector)
 
template<typename Number2 >
Vector< Number > & operator= (const Vector< Number2 > &in_vector)
 
Vector< Number > & operator= (const Number s)
 
virtual Vector< Number > & operator*= (const Number factor)
 
virtual Vector< Number > & operator/= (const Number factor)
 
virtual Vector< Number > & operator+= (const VectorSpaceVector< Number > &V)
 
virtual Vector< Number > & operator-= (const VectorSpaceVector< Number > &V)
 
virtual Number operator* (const VectorSpaceVector< Number > &V) const
 
virtual void import (const ReadWriteVector< Number > &V, VectorOperation::values operation, std_cxx11::shared_ptr< const CommunicationPatternBase > communication_pattern=std_cxx11::shared_ptr< const CommunicationPatternBase >())
 
virtual void add (const Number a)
 
virtual void add (const Number a, const VectorSpaceVector< Number > &V)
 
virtual void add (const Number a, const VectorSpaceVector< Number > &V, const Number b, const VectorSpaceVector< Number > &W)
 
virtual void sadd (const Number s, const Number a, const VectorSpaceVector< Number > &V)
 
virtual void scale (const VectorSpaceVector< Number > &scaling_factors)
 
virtual void equ (const Number a, const VectorSpaceVector< Number > &V)
 
virtual value_type mean_value () const
 
virtual VectorSpaceVector< Number >::real_type l1_norm () const
 
virtual VectorSpaceVector< Number >::real_type l2_norm () const
 
virtual VectorSpaceVector< Number >::real_type linfty_norm () const
 
virtual Number add_and_dot (const Number a, const VectorSpaceVector< Number > &V, const VectorSpaceVector< Number > &W)
 
virtual size_type size () const
 
virtual ::IndexSet locally_owned_elements () const
 
virtual void print (std::ostream &out, const unsigned int precision=3, const bool scientific=true, const bool across=true) const
 
void block_write (std::ostream &out) const
 
void block_read (std::istream &in)
 
virtual std::size_t memory_consumption () const
 
- 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 ()
 
void reinit (const size_type size, const bool omit_zeroing_entries=false)
 
template<typename Number2 >
void reinit (const ReadWriteVector< Number2 > &in_vector, const bool omit_zeroing_entries=false)
 
void reinit (const IndexSet &locally_stored_indices, const bool omit_zeroing_entries=false)
 
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, std_cxx11::shared_ptr< const CommunicationPatternBase > communication_pattern=std_cxx11::shared_ptr< const CommunicationPatternBase >())
 
void import (const PETScWrappers::MPI::Vector &petsc_vec, VectorOperation::values operation, std_cxx11::shared_ptr< const CommunicationPatternBase > communication_pattern=std_cxx11::shared_ptr< const CommunicationPatternBase >())
 
void import (const TrilinosWrappers::MPI::Vector &trilinos_vec, VectorOperation::values operation, std_cxx11::shared_ptr< const CommunicationPatternBase > communication_pattern=std_cxx11::shared_ptr< const CommunicationPatternBase >())
 
void import (const EpetraWrappers::Vector &epetra_vec, VectorOperation::values operation, std_cxx11::shared_ptr< const CommunicationPatternBase > communication_pattern=std_cxx11::shared_ptr< const CommunicationPatternBase >())
 
size_type size () const
 
size_type n_elements () const
 
const IndexSetget_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 &&)
 
virtual ~Subscriptor ()
 
Subscriptoroperator= (const Subscriptor &)
 
Subscriptoroperator= (Subscriptor &&)
 
void subscribe (const char *identifier=0) const
 
void unsubscribe (const char *identifier=0) 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 ()
 

Static Public Member Functions

static ::ExceptionBaseExcVectorTypeNotCompatible ()
 
- Static Public Member Functions inherited from Subscriptor
static ::ExceptionBaseExcInUse (int arg1, char *arg2, std::string &arg3)
 
static ::ExceptionBaseExcNoSubscriber (char *arg1, char *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, std_cxx11::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_cxx11::shared_ptr< CommunicationPatternBasecomm_pattern
 
Number * val
 
std_cxx11::shared_ptr< ::parallel::internal::TBBPartitioner > thread_loop_partitioner
 

Detailed Description

template<typename Number>
class LinearAlgebra::Vector< Number >

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.

Author
Bruno Turcksin, 2015.

Definition at line 62 of file la_vector.h.

Constructor & Destructor Documentation

◆ Vector() [1/4]

template<typename Number >
Vector< Number >::Vector ( )
inline

Constructor. Create a vector of dimension zero.

Definition at line 301 of file la_vector.h.

◆ Vector() [2/4]

template<typename Number >
Vector< Number >::Vector ( const Vector< Number > &  V)
inline

Copy constructor. Sets the dimension to that of the given vector and copies all elements.

Definition at line 307 of file la_vector.h.

◆ Vector() [3/4]

template<typename Number >
Vector< Number >::Vector ( const size_type  n)
inlineexplicit

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 316 of file la_vector.h.

◆ Vector() [4/4]

template<typename Number >
template<typename InputIterator >
Vector< Number >::Vector ( const InputIterator  first,
const InputIterator  last 
)
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 326 of file la_vector.h.

Member Function Documentation

◆ operator=() [1/3]

template<typename Number>
Vector<Number>& LinearAlgebra::Vector< Number >::operator= ( const Vector< Number > &  in_vector)

Copies the data of the input vector in_vector.

◆ operator=() [2/3]

template<typename Number>
template<typename Number2 >
Vector<Number>& LinearAlgebra::Vector< Number >::operator= ( const Vector< Number2 > &  in_vector)

Copies the data of the input vector in_vector.

◆ operator=() [3/3]

template<typename Number>
Vector<Number>& LinearAlgebra::Vector< Number >::operator= ( const Number  s)

Sets all elements of the vector to the scalar s. This operation is only allowed if s is equal to zero.

◆ operator*=()

template<typename Number>
virtual Vector<Number>& LinearAlgebra::Vector< Number >::operator*= ( const Number  factor)
virtual

Multiply the entire vector by a fixed factor.

Implements LinearAlgebra::VectorSpaceVector< Number >.

◆ operator/=()

template<typename Number>
virtual Vector<Number>& LinearAlgebra::Vector< Number >::operator/= ( const Number  factor)
virtual

Divide the entire vector by a fixed factor.

Implements LinearAlgebra::VectorSpaceVector< Number >.

◆ operator+=()

template<typename Number>
virtual Vector<Number>& LinearAlgebra::Vector< Number >::operator+= ( const VectorSpaceVector< Number > &  V)
virtual

Add the vector V to the present one.

Implements LinearAlgebra::VectorSpaceVector< Number >.

◆ operator-=()

template<typename Number>
virtual Vector<Number>& LinearAlgebra::Vector< Number >::operator-= ( const VectorSpaceVector< Number > &  V)
virtual

Subtract the vector V from the present one.

Implements LinearAlgebra::VectorSpaceVector< Number >.

◆ operator*()

template<typename Number>
virtual Number LinearAlgebra::Vector< Number >::operator* ( const VectorSpaceVector< Number > &  V) const
virtual

Return the scalar product of two vectors.

Implements LinearAlgebra::VectorSpaceVector< Number >.

◆ import()

template<typename Number>
virtual void LinearAlgebra::Vector< Number >::import ( const ReadWriteVector< Number > &  V,
VectorOperation::values  operation,
std_cxx11::shared_ptr< const CommunicationPatternBase communication_pattern = std_cxx11::shared_ptr< const CommunicationPatternBase >() 
)
virtual

This function is not implemented and will throw an exception.

Implements LinearAlgebra::VectorSpaceVector< Number >.

◆ add() [1/3]

template<typename Number>
virtual void LinearAlgebra::Vector< Number >::add ( const Number  a)
virtual

Add a to all components. Note that a is a scalar not a vector.

Implements LinearAlgebra::VectorSpaceVector< Number >.

◆ add() [2/3]

template<typename Number>
virtual void LinearAlgebra::Vector< Number >::add ( const Number  a,
const VectorSpaceVector< Number > &  V 
)
virtual

Simple addition of a multiple of a vector, i.e. *this += a*V.

Implements LinearAlgebra::VectorSpaceVector< Number >.

◆ add() [3/3]

template<typename Number>
virtual void LinearAlgebra::Vector< Number >::add ( const Number  a,
const VectorSpaceVector< Number > &  V,
const Number  b,
const VectorSpaceVector< Number > &  W 
)
virtual

Multiple addition of a multiple of a vector, i.e. *this += a*V+b*W.

Implements LinearAlgebra::VectorSpaceVector< Number >.

◆ sadd()

template<typename Number>
virtual void LinearAlgebra::Vector< Number >::sadd ( const Number  s,
const Number  a,
const VectorSpaceVector< Number > &  V 
)
virtual

Scaling and simple addition of a multiple of a vector, i.e. *this = s*(*this)+a*V.

Implements LinearAlgebra::VectorSpaceVector< Number >.

◆ scale()

template<typename Number>
virtual void LinearAlgebra::Vector< Number >::scale ( const VectorSpaceVector< Number > &  scaling_factors)
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.

Implements LinearAlgebra::VectorSpaceVector< Number >.

◆ equ()

template<typename Number>
virtual void LinearAlgebra::Vector< Number >::equ ( const Number  a,
const VectorSpaceVector< Number > &  V 
)
virtual

Assignment *this = a*V.

Implements LinearAlgebra::VectorSpaceVector< Number >.

◆ mean_value()

template<typename Number>
virtual value_type LinearAlgebra::Vector< Number >::mean_value ( ) const
virtual

Return the mean value of all the entries of this vector.

Implements LinearAlgebra::VectorSpaceVector< Number >.

◆ l1_norm()

template<typename Number>
virtual VectorSpaceVector<Number>::real_type LinearAlgebra::Vector< Number >::l1_norm ( ) const
virtual

Return the l1 norm of the vector (i.e., the sum of the absolute values of all entries).

Implements LinearAlgebra::VectorSpaceVector< Number >.

◆ l2_norm()

template<typename Number>
virtual VectorSpaceVector<Number>::real_type LinearAlgebra::Vector< Number >::l2_norm ( ) const
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).

Implements LinearAlgebra::VectorSpaceVector< Number >.

◆ linfty_norm()

template<typename Number>
virtual VectorSpaceVector<Number>::real_type LinearAlgebra::Vector< Number >::linfty_norm ( ) const
virtual

Return the maximum norm of the vector (i.e., the maximum absolute value among all entries and among all processors).

Implements LinearAlgebra::VectorSpaceVector< Number >.

◆ add_and_dot()

template<typename Number>
virtual Number LinearAlgebra::Vector< Number >::add_and_dot ( const Number  a,
const VectorSpaceVector< Number > &  V,
const VectorSpaceVector< Number > &  W 
)
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

this->add(a, V);
return_value = *this * W;

Implements LinearAlgebra::VectorSpaceVector< Number >.

◆ size()

template<typename Number >
Vector< Number >::size_type Vector< Number >::size ( ) const
inlinevirtual

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 336 of file la_vector.h.

◆ locally_owned_elements()

template<typename Number >
IndexSet Vector< Number >::locally_owned_elements ( ) const
inlinevirtual

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 procesors 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

vec.locally_owned_elements() == complete_index_set(vec.size())

Implements LinearAlgebra::VectorSpaceVector< Number >.

Definition at line 345 of file la_vector.h.

◆ print()

template<typename Number >
void Vector< Number >::print ( std::ostream &  out,
const unsigned int  precision = 3,
const bool  scientific = true,
const bool  across = true 
) const
inlinevirtual

Prints the vector to the output stream out.

Implements LinearAlgebra::VectorSpaceVector< Number >.

Definition at line 354 of file la_vector.h.

◆ block_write()

template<typename Number>
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.

◆ block_read()

template<typename Number>
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.

◆ memory_consumption()

template<typename Number >
std::size_t Vector< Number >::memory_consumption ( ) const
inlinevirtual

Return the memory consumption of this class in bytes.

Implements LinearAlgebra::VectorSpaceVector< Number >.

Definition at line 382 of file la_vector.h.

◆ serialize()

template<typename Number >
template<typename Archive >
void Vector< Number >::serialize ( Archive &  ar,
const unsigned int  version 
)
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 367 of file la_vector.h.

Friends And Related Function Documentation

◆ Vector

template<typename Number>
template<typename Number2 >
friend class Vector
friend

Make all other ReadWriteVector types friends.

Typedef for the vector type used.

Definition at line 293 of file la_vector.h.


The documentation for this class was generated from the following file: