Reference documentation for deal.II version 8.5.1
Public Types | Public Member Functions | Static Public Attributes | Protected Member Functions | Related Functions | List of all members

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

Inheritance diagram for PETScWrappers::Vector:
[legend]

Public Types

typedef types::global_dof_index size_type
 
- Public Types inherited from PETScWrappers::VectorBase
typedef PetscScalar value_type
 

Public Member Functions

 Vector ()
 
 Vector (const size_type n)
 
template<typename Number >
 Vector (const ::Vector< Number > &v)
 
 Vector (const Vec &v)
 
 Vector (const Vector &v)
 
 Vector (const MPI::Vector &v)
 
void clear ()
 
Vectoroperator= (const Vector &v)
 
Vectoroperator= (const MPI::Vector &v)
 
Vectoroperator= (const PetscScalar s)
 
template<typename number >
Vectoroperator= (const ::Vector< number > &v)
 
void reinit (const size_type N, const bool omit_zeroing_entries=false)
 
void reinit (const Vector &v, const bool omit_zeroing_entries=false)
 
- Public Member Functions inherited from PETScWrappers::VectorBase
 VectorBase ()
 
 VectorBase (const VectorBase &v)
 
 VectorBase (const Vec &v)
 
virtual ~VectorBase ()
 
void compress (const VectorOperation::values operation)
 
VectorBaseoperator= (const PetscScalar s)
 
bool operator== (const VectorBase &v) const
 
bool operator!= (const VectorBase &v) const
 
size_type size () const
 
size_type local_size () const
 
std::pair< size_type, size_type > local_range () const
 
bool in_local_range (const size_type index) const
 
IndexSet locally_owned_elements () const
 
bool has_ghost_elements () const
 
reference operator() (const size_type index)
 
PetscScalar operator() (const size_type index) const
 
reference operator[] (const size_type index)
 
PetscScalar operator[] (const size_type index) const
 
void set (const std::vector< size_type > &indices, const std::vector< PetscScalar > &values)
 
void extract_subvector_to (const std::vector< size_type > &indices, std::vector< PetscScalar > &values) const
 
template<typename ForwardIterator , typename OutputIterator >
void extract_subvector_to (const ForwardIterator indices_begin, const ForwardIterator indices_end, OutputIterator values_begin) const
 
void add (const std::vector< size_type > &indices, const std::vector< PetscScalar > &values)
 
void add (const std::vector< size_type > &indices, const ::Vector< PetscScalar > &values)
 
void add (const size_type n_elements, const size_type *indices, const PetscScalar *values)
 
PetscScalar operator* (const VectorBase &vec) const
 
real_type norm_sqr () const
 
PetscScalar mean_value () const
 
real_type l1_norm () const
 
real_type l2_norm () const
 
real_type lp_norm (const real_type p) const
 
real_type linfty_norm () const
 
PetscScalar add_and_dot (const PetscScalar a, const VectorBase &V, const VectorBase &W)
 
real_type normalize () const 1
 
real_type min () const
 
real_type max () const
 
VectorBaseabs () 1
 
VectorBaseconjugate () 1
 
VectorBasemult () 1
 
VectorBasemult (const VectorBase &v) 1
 
VectorBasemult (const VectorBase &u, const VectorBase &v) 1
 
bool all_zero () const
 
bool is_non_negative () const
 
VectorBaseoperator*= (const PetscScalar factor)
 
VectorBaseoperator/= (const PetscScalar factor)
 
VectorBaseoperator+= (const VectorBase &V)
 
VectorBaseoperator-= (const VectorBase &V)
 
void add (const PetscScalar s)
 
void add (const VectorBase &V) 1
 
void add (const PetscScalar a, const VectorBase &V)
 
void add (const PetscScalar a, const VectorBase &V, const PetscScalar b, const VectorBase &W)
 
void sadd (const PetscScalar s, const VectorBase &V)
 
void sadd (const PetscScalar s, const PetscScalar a, const VectorBase &V)
 
void sadd (const PetscScalar s, const PetscScalar a, const VectorBase &V, const PetscScalar b, const VectorBase &W) 1
 
void sadd (const PetscScalar s, const PetscScalar a, const VectorBase &V, const PetscScalar b, const VectorBase &W, const PetscScalar c, const VectorBase &X) 1
 
void scale (const VectorBase &scaling_factors)
 
void equ (const PetscScalar a, const VectorBase &V)
 
void equ (const PetscScalar a, const VectorBase &V, const PetscScalar b, const VectorBase &W) 1
 
void ratio (const VectorBase &a, const VectorBase &b) 1
 
void write_ascii (const PetscViewerFormat format=PETSC_VIEWER_DEFAULT)
 
void print (std::ostream &out, const unsigned int precision=3, const bool scientific=true, const bool across=true) const
 
void swap (VectorBase &v)
 
 operator const Vec & () const
 
std::size_t memory_consumption () const
 
virtual const MPI_Comm & get_mpi_communicator () 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)
 

Static Public Attributes

static const bool supports_distributed_data = false
 

Protected Member Functions

void create_vector (const size_type n)
 
- Protected Member Functions inherited from PETScWrappers::VectorBase
void do_set_add_operation (const size_type n_elements, const size_type *indices, const PetscScalar *values, const bool add_values)
 

Related Functions

(Note that these are not member functions.)

void swap (Vector &u, Vector &v)
 

Additional Inherited Members

- Static Public Member Functions inherited from Subscriptor
static ::ExceptionBaseExcInUse (int arg1, char *arg2, std::string &arg3)
 
static ::ExceptionBaseExcNoSubscriber (char *arg1, char *arg2)
 
- Protected Attributes inherited from PETScWrappers::VectorBase
Vec vector
 
bool ghosted
 
IndexSet ghost_indices
 
VectorOperation::values last_action
 
bool attained_ownership
 

Detailed Description

Implementation of a sequential vector class based on PETSC. All the functionality is actually in the base class, except for the calls to generate a sequential vector. This is possible since PETSc only works on an abstract vector type and internally distributes to functions that do the actual work depending on the actual vector type (much like using virtual functions). Only the functions creating a vector of specific type differ, and are implemented in this particular class.

This class is deprecated, use PETScWrappers::MPI::Vector instead.

Author
Wolfgang Bangerth, 2004

Definition at line 54 of file petsc_vector.h.

Member Typedef Documentation

◆ size_type

Declare type for container size.

Definition at line 60 of file petsc_vector.h.

Constructor & Destructor Documentation

◆ Vector() [1/6]

Vector< Number >::Vector ( )

Default constructor. Initialize the vector as empty.

Definition at line 30 of file petsc_vector.cc.

◆ Vector() [2/6]

Vector< Number >::Vector ( const size_type  n)
explicit

Constructor. Set dimension to n and initialize all elements with zero.

The constructor is made explicit to avoid accidents 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 37 of file petsc_vector.cc.

◆ Vector() [3/6]

template<typename Number >
PETScWrappers::Vector::Vector ( const ::Vector< Number > &  v)
explicit

Copy-constructor from deal.II vectors. Sets the dimension to that of the given vector, and copies all elements.

◆ Vector() [4/6]

PETScWrappers::Vector::Vector ( const Vec &  v)
explicit

Construct it from an existing PETSc Vector of type Vec. Note: this does not copy the contents and just keeps a pointer. You need to make sure the vector is not used twice at the same time or destroyed while in use. This class does not destroy the PETSc object. Handle with care!

◆ Vector() [5/6]

Vector< Number >::Vector ( const Vector v)

Copy-constructor the values from a PETSc wrapper vector class.

Definition at line 44 of file petsc_vector.cc.

◆ Vector() [6/6]

Vector< Number >::Vector ( const MPI::Vector v)
explicit

Copy-constructor: copy the values from a PETSc wrapper parallel vector class.

Note that due to the communication model of MPI, all processes have to actually perform this operation, even if they do not use the result. It is not sufficient if only one processor tries to copy the elements from the other processors over to its own process space.

Definition at line 56 of file petsc_vector.cc.

Member Function Documentation

◆ clear()

void Vector< Number >::clear ( )
virtual

Release all memory and return to a state just like after having called the default constructor.

Reimplemented from PETScWrappers::VectorBase.

Definition at line 67 of file petsc_vector.cc.

◆ operator=() [1/4]

Vector& PETScWrappers::Vector::operator= ( const Vector v)

Copy the given vector. Resize the present vector if necessary.

◆ operator=() [2/4]

Vector& PETScWrappers::Vector::operator= ( const MPI::Vector v)

Copy all the elements of the parallel vector v into this local vector. Note that due to the communication model of MPI, all processes have to actually perform this operation, even if they do not use the result. It is not sufficient if only one processor tries to copy the elements from the other processors over to its own process space.

◆ operator=() [3/4]

Vector& PETScWrappers::Vector::operator= ( const PetscScalar  s)

Set all components of the vector to the given number s. Simply pass this down to the base class, but we still need to declare this function to make the example given in the discussion about making the constructor explicit work.

Since the semantics of assigning a scalar to a vector are not immediately clear, this operator should really only be used if you want to set the entire vector to zero. This allows the intuitive notation v=0. Assigning other values is deprecated and may be disallowed in the future.

◆ operator=() [4/4]

template<typename number >
Vector& PETScWrappers::Vector::operator= ( const ::Vector< number > &  v)

Copy the values of a deal.II vector (as opposed to those of the PETSc vector wrapper class) into this object.

◆ reinit() [1/2]

void Vector< Number >::reinit ( const size_type  N,
const bool  omit_zeroing_entries = false 
)

Change the dimension of the vector to N. It is unspecified how resizing the vector affects the memory allocation of this object; i.e., it is not guaranteed that resizing it to a smaller size actually also reduces memory consumption, or if for efficiency the same amount of memory is used for less data.

If omit_zeroing_entries is false, the vector is filled by zeros. Otherwise, the elements are left an unspecified state.

Definition at line 76 of file petsc_vector.cc.

◆ reinit() [2/2]

void Vector< Number >::reinit ( const Vector v,
const bool  omit_zeroing_entries = false 
)

Change the dimension to that of the vector v. The same applies as for the other reinit() function.

The elements of v are not copied, i.e. this function is the same as calling reinit (v.size(), omit_zeroing_entries).

Definition at line 113 of file petsc_vector.cc.

◆ create_vector()

void Vector< Number >::create_vector ( const size_type  n)
protected

Create a vector of length n. For this class, we create a sequential vector. n denotes the total size of the vector to be created.

Definition at line 122 of file petsc_vector.cc.

Friends And Related Function Documentation

◆ swap()

void swap ( Vector u,
Vector v 
)
related

Global function swap which overloads the default implementation of the C++ standard library which uses a temporary object. The function simply exchanges the data of the two vectors.

Author
Wolfgang Bangerth, 2004

Definition at line 212 of file petsc_vector.h.

Member Data Documentation

◆ supports_distributed_data

const bool PETScWrappers::Vector::supports_distributed_data = false
static

A variable that indicates whether this vector supports distributed data storage. If true, then this vector also needs an appropriate compress() function that allows communicating recent set or add operations to individual elements to be communicated to other processors.

For the current class, the variable equals false, since it does not support parallel data storage. If you do need parallel data storage, use PETScWrappers::MPI::Vector.

Deprecated:
instead of using this variable, please use the type trait value is_serial_vector< VectorType >::value

Definition at line 75 of file petsc_vector.h.


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