Reference documentation for deal.II version 8.5.1
Public Types | Public Member Functions | Protected Member Functions | Protected Attributes | Private Member Functions | Friends | Related Functions | List of all members
PETScWrappers::VectorBase Class Reference

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

Inheritance diagram for PETScWrappers::VectorBase:
[legend]

Public Types

typedef PetscScalar value_type
 

Public Member Functions

 VectorBase ()
 
 VectorBase (const VectorBase &v)
 
 VectorBase (const Vec &v)
 
virtual ~VectorBase ()
 
virtual void clear ()
 
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)
 

Protected Member Functions

void do_set_add_operation (const size_type n_elements, const size_type *indices, const PetscScalar *values, const bool add_values)
 

Protected Attributes

Vec vector
 
bool ghosted
 
IndexSet ghost_indices
 
VectorOperation::values last_action
 
bool attained_ownership
 

Private Member Functions

VectorBaseoperator= (const VectorBase &)
 

Friends

class internal::VectorReference
 

Related Functions

(Note that these are not member functions.)

void swap (VectorBase &u, VectorBase &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)
 

Detailed Description

Base class for all vector classes that are implemented on top of the PETSc vector types. Since in PETSc all vector types (i.e. sequential and parallel ones) are built by filling the contents of an abstract object that is only referenced through a pointer of a type that is independent of the actual vector type, we can implement almost all functionality of vectors in this base class. Derived classes will then only have to provide the functionality to create one or the other kind of vector.

The interface of this class is modeled after the existing Vector class in deal.II. It has almost the same member functions, and is often exchangeable. However, since PETSc only supports a single scalar type (either double, float, or a complex data type), it is not templated, and only works with whatever your PETSc installation has defined the data type PetscScalar to.

Note that PETSc only guarantees that operations do what you expect if the functions VecAssemblyBegin and VecAssemblyEnd have been called after vector assembly. Therefore, you need to call Vector::compress() before you actually use the vector.

Author
Wolfgang Bangerth, 2004

Definition at line 230 of file petsc_vector_base.h.

Member Typedef Documentation

◆ value_type

Declare some of the standard types used in all containers. These types parallel those in the C++ standard libraries vector<...> class.

Definition at line 238 of file petsc_vector_base.h.

Constructor & Destructor Documentation

◆ VectorBase() [1/3]

PETScWrappers::VectorBase::VectorBase ( )

Default constructor. It doesn't do anything, derived classes will have to initialize the data.

Definition at line 150 of file petsc_vector_base.cc.

◆ VectorBase() [2/3]

PETScWrappers::VectorBase::VectorBase ( const VectorBase v)

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

Definition at line 164 of file petsc_vector_base.cc.

◆ VectorBase() [3/3]

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

Initialize a Vector from a PETSc Vec object. Note that we do not copy the vector and we do not attain ownership, so we do not destroy the PETSc object in the destructor.

Definition at line 185 of file petsc_vector_base.cc.

◆ ~VectorBase()

PETScWrappers::VectorBase::~VectorBase ( )
virtual

Destructor

Definition at line 200 of file petsc_vector_base.cc.

Member Function Documentation

◆ clear()

void PETScWrappers::VectorBase::clear ( )
virtual

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

Reimplemented in PETScWrappers::MPI::Vector, and PETScWrappers::Vector.

Definition at line 217 of file petsc_vector_base.cc.

◆ compress()

void PETScWrappers::VectorBase::compress ( const VectorOperation::values  operation)

Compress the underlying representation of the PETSc object, i.e. flush the buffers of the vector object if it has any. This function is necessary after writing into a vector element-by-element and before anything else can be done on it.

See Compressing distributed objects for more information.

Definition at line 407 of file petsc_vector_base.cc.

◆ operator=() [1/2]

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

Set all components of the vector to the given number s. Simply pass this down to the individual block objects, 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.

Definition at line 238 of file petsc_vector_base.cc.

◆ operator==()

bool PETScWrappers::VectorBase::operator== ( const VectorBase v) const

Test for equality. This function assumes that the present vector and the one to compare with have the same size already, since comparing vectors of different sizes makes not much sense anyway.

Definition at line 266 of file petsc_vector_base.cc.

◆ operator!=()

bool PETScWrappers::VectorBase::operator!= ( const VectorBase v) const

Test for inequality. This function assumes that the present vector and the one to compare with have the same size already, since comparing vectors of different sizes makes not much sense anyway.

Definition at line 281 of file petsc_vector_base.cc.

◆ size()

VectorBase::size_type PETScWrappers::VectorBase::size ( ) const

Return the global dimension of the vector.

Definition at line 296 of file petsc_vector_base.cc.

◆ local_size()

VectorBase::size_type PETScWrappers::VectorBase::local_size ( ) const

Return the local dimension of the vector, i.e. the number of elements stored on the present MPI process. For sequential vectors, this number is the same as size(), but for parallel vectors it may be smaller.

To figure out which elements exactly are stored locally, use local_range().

Definition at line 308 of file petsc_vector_base.cc.

◆ local_range()

std::pair< VectorBase::size_type, VectorBase::size_type > PETScWrappers::VectorBase::local_range ( ) const

Return a pair of indices indicating which elements of this vector are stored locally. The first number is the index of the first element stored, the second the index of the one past the last one that is stored locally. If this is a sequential vector, then the result will be the pair (0,N), otherwise it will be a pair (i,i+n), where n=local_size().

Definition at line 320 of file petsc_vector_base.cc.

◆ in_local_range()

bool PETScWrappers::VectorBase::in_local_range ( const size_type  index) const

Return whether index is in the local range or not, see also local_range().

◆ locally_owned_elements()

IndexSet PETScWrappers::VectorBase::locally_owned_elements ( ) const

Return an index set that describes which elements of this vector are owned by the current processor. Note that this index set does not include elements this vector may store locally as ghost elements but that are in fact owned by another 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

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

◆ has_ghost_elements()

bool PETScWrappers::VectorBase::has_ghost_elements ( ) const

Return if the vector contains ghost elements.

See also
vectors with ghost elements

◆ operator()() [1/2]

reference PETScWrappers::VectorBase::operator() ( const size_type  index)

Provide access to a given element, both read and write.

◆ operator()() [2/2]

PetscScalar PETScWrappers::VectorBase::operator() ( const size_type  index) const

Provide read-only access to an element.

◆ operator[]() [1/2]

reference PETScWrappers::VectorBase::operator[] ( const size_type  index)

Provide access to a given element, both read and write.

Exactly the same as operator().

◆ operator[]() [2/2]

PetscScalar PETScWrappers::VectorBase::operator[] ( const size_type  index) const

Provide read-only access to an element.

Exactly the same as operator().

◆ set()

void PETScWrappers::VectorBase::set ( const std::vector< size_type > &  indices,
const std::vector< PetscScalar > &  values 
)

A collective set operation: instead of setting individual elements of a vector, this function allows to set a whole set of elements at once. The indices of the elements to be set are stated in the first argument, the corresponding values in the second.

Definition at line 333 of file petsc_vector_base.cc.

◆ extract_subvector_to() [1/2]

void PETScWrappers::VectorBase::extract_subvector_to ( const std::vector< size_type > &  indices,
std::vector< PetscScalar > &  values 
) const

A collective get operation: instead of getting individual elements of a vector, this function allows to get a whole set of elements at once. The indices of the elements to be read are stated in the first argument, the corresponding values are returned in the second.

◆ extract_subvector_to() [2/2]

template<typename ForwardIterator , typename OutputIterator >
void PETScWrappers::VectorBase::extract_subvector_to ( const ForwardIterator  indices_begin,
const ForwardIterator  indices_end,
OutputIterator  values_begin 
) const

Just as the above, but with pointers. Useful in minimizing copying of data around.

◆ add() [1/7]

void PETScWrappers::VectorBase::add ( const std::vector< size_type > &  indices,
const std::vector< PetscScalar > &  values 
)

A collective add operation: This function adds a whole set of values stored in values to the vector components specified by indices.

Definition at line 344 of file petsc_vector_base.cc.

◆ add() [2/7]

void PETScWrappers::VectorBase::add ( const std::vector< size_type > &  indices,
const ::Vector< PetscScalar > &  values 
)

This is a second collective add operation. As a difference, this function takes a deal.II vector of values.

Definition at line 355 of file petsc_vector_base.cc.

◆ add() [3/7]

void PETScWrappers::VectorBase::add ( const size_type  n_elements,
const size_type *  indices,
const PetscScalar *  values 
)

Take an address where n_elements are stored contiguously and add them into the vector. Handles all cases which are not covered by the other two add() functions above.

Definition at line 366 of file petsc_vector_base.cc.

◆ operator*()

PetscScalar PETScWrappers::VectorBase::operator* ( const VectorBase vec) const

Return the scalar product of two vectors. The vectors must have the same size.

For complex valued vector, this gives \(\left(v^\ast,vec\right)\).

Definition at line 376 of file petsc_vector_base.cc.

◆ norm_sqr()

VectorBase::real_type PETScWrappers::VectorBase::norm_sqr ( ) const

Return square of the \(l_2\)-norm.

Definition at line 462 of file petsc_vector_base.cc.

◆ mean_value()

PetscScalar PETScWrappers::VectorBase::mean_value ( ) const

Return the mean value of the elements of this vector.

Definition at line 471 of file petsc_vector_base.cc.

◆ l1_norm()

VectorBase::real_type PETScWrappers::VectorBase::l1_norm ( ) const

\(l_1\)-norm of the vector. The sum of the absolute values.

Note
In complex-valued PETSc priori to 3.7.0 this norm is implemented as the sum of absolute values of real and imaginary parts of elements of a complex vector.

Definition at line 525 of file petsc_vector_base.cc.

◆ l2_norm()

VectorBase::real_type PETScWrappers::VectorBase::l2_norm ( ) const

\(l_2\)-norm of the vector. The square root of the sum of the squares of the elements.

Definition at line 538 of file petsc_vector_base.cc.

◆ lp_norm()

VectorBase::real_type PETScWrappers::VectorBase::lp_norm ( const real_type  p) const

\(l_p\)-norm of the vector. The pth root of the sum of the pth powers of the absolute values of the elements.

Definition at line 551 of file petsc_vector_base.cc.

◆ linfty_norm()

VectorBase::real_type PETScWrappers::VectorBase::linfty_norm ( ) const

\(l_\infty\)-norm of the vector. Return the value of the vector element with the maximum absolute value.

Definition at line 597 of file petsc_vector_base.cc.

◆ add_and_dot()

PetscScalar PETScWrappers::VectorBase::add_and_dot ( const PetscScalar  a,
const VectorBase V,
const VectorBase W 
)

Performs 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;

The reason this function exists is for compatibility with deal.II's own vector classes which can implement this functionality with less memory transfer. However, for PETSc vectors such a combined operation is not natively supported and thus the cost is completely equivalent as calling the two methods separately.

Definition at line 396 of file petsc_vector_base.cc.

◆ normalize()

VectorBase::real_type PETScWrappers::VectorBase::normalize ( ) const

Normalize vector by dividing by the \(l_2\)-norm of the vector. Return the vector norm before normalization.

This function is deprecated.

Definition at line 610 of file petsc_vector_base.cc.

◆ min()

VectorBase::real_type PETScWrappers::VectorBase::min ( ) const

Return the value of the vector element with the largest negative value.

Definition at line 622 of file petsc_vector_base.cc.

◆ max()

VectorBase::real_type PETScWrappers::VectorBase::max ( ) const

Return the value of the vector element with the largest positive value.

Definition at line 635 of file petsc_vector_base.cc.

◆ abs()

VectorBase & PETScWrappers::VectorBase::abs ( )

Replace every element in a vector with its absolute value.

This function is deprecated.

Definition at line 648 of file petsc_vector_base.cc.

◆ conjugate()

VectorBase & PETScWrappers::VectorBase::conjugate ( )

Conjugate a vector.

This function is deprecated.

Definition at line 661 of file petsc_vector_base.cc.

◆ mult() [1/3]

VectorBase & PETScWrappers::VectorBase::mult ( )

A collective piecewise multiply operation on this vector with itself. TODO: The model for this function should be similar to add ().

This function is deprecated.

Definition at line 674 of file petsc_vector_base.cc.

◆ mult() [2/3]

VectorBase & PETScWrappers::VectorBase::mult ( const VectorBase v)

Same as above, but a collective piecewise multiply operation of this vector with v.

This function is deprecated.

Definition at line 686 of file petsc_vector_base.cc.

◆ mult() [3/3]

VectorBase & PETScWrappers::VectorBase::mult ( const VectorBase u,
const VectorBase v 
)

Same as above, but a collective piecewise multiply operation of u with v.

This function is deprecated.

Definition at line 697 of file petsc_vector_base.cc.

◆ all_zero()

bool PETScWrappers::VectorBase::all_zero ( ) const

Return whether the vector contains only elements with value zero. This is a collective operation. This function is expensive, because potentially all elements have to be checked.

Definition at line 708 of file petsc_vector_base.cc.

◆ is_non_negative()

bool PETScWrappers::VectorBase::is_non_negative ( ) const

Return true if the vector has no negative entries, i.e. all entries are zero or positive. This function is used, for example, to check whether refinement indicators are really all positive (or zero).

Definition at line 761 of file petsc_vector_base.cc.

◆ operator*=()

VectorBase & PETScWrappers::VectorBase::operator*= ( const PetscScalar  factor)

Multiply the entire vector by a fixed factor.

Definition at line 793 of file petsc_vector_base.cc.

◆ operator/=()

VectorBase & PETScWrappers::VectorBase::operator/= ( const PetscScalar  factor)

Divide the entire vector by a fixed factor.

Definition at line 807 of file petsc_vector_base.cc.

◆ operator+=()

VectorBase & PETScWrappers::VectorBase::operator+= ( const VectorBase V)

Add the given vector to the present one.

Definition at line 824 of file petsc_vector_base.cc.

◆ operator-=()

VectorBase & PETScWrappers::VectorBase::operator-= ( const VectorBase V)

Subtract the given vector from the present one.

Definition at line 836 of file petsc_vector_base.cc.

◆ add() [4/7]

void PETScWrappers::VectorBase::add ( const PetscScalar  s)

Addition of s to all components. Note that s is a scalar and not a vector.

Definition at line 848 of file petsc_vector_base.cc.

◆ add() [5/7]

void PETScWrappers::VectorBase::add ( const VectorBase V)

Simple vector addition, equal to the operator +=.

Deprecated:
Use the operator += instead.

Definition at line 860 of file petsc_vector_base.cc.

◆ add() [6/7]

void PETScWrappers::VectorBase::add ( const PetscScalar  a,
const VectorBase V 
)

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

Definition at line 868 of file petsc_vector_base.cc.

◆ add() [7/7]

void PETScWrappers::VectorBase::add ( const PetscScalar  a,
const VectorBase V,
const PetscScalar  b,
const VectorBase W 
)

Multiple addition of scaled vectors, i.e. *this += a*V+b*W.

Definition at line 881 of file petsc_vector_base.cc.

◆ sadd() [1/4]

void PETScWrappers::VectorBase::sadd ( const PetscScalar  s,
const VectorBase V 
)

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

Definition at line 900 of file petsc_vector_base.cc.

◆ sadd() [2/4]

void PETScWrappers::VectorBase::sadd ( const PetscScalar  s,
const PetscScalar  a,
const VectorBase V 
)

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

Definition at line 913 of file petsc_vector_base.cc.

◆ sadd() [3/4]

void PETScWrappers::VectorBase::sadd ( const PetscScalar  s,
const PetscScalar  a,
const VectorBase V,
const PetscScalar  b,
const VectorBase W 
)

Scaling and multiple addition.

This function is deprecated.

Definition at line 931 of file petsc_vector_base.cc.

◆ sadd() [4/4]

void PETScWrappers::VectorBase::sadd ( const PetscScalar  s,
const PetscScalar  a,
const VectorBase V,
const PetscScalar  b,
const VectorBase W,
const PetscScalar  c,
const VectorBase X 
)

Scaling and multiple addition. *this = s*(*this)+a*V + b*W + c*X.

This function is deprecated.

Definition at line 956 of file petsc_vector_base.cc.

◆ scale()

void PETScWrappers::VectorBase::scale ( const VectorBase scaling_factors)

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 984 of file petsc_vector_base.cc.

◆ equ() [1/2]

void PETScWrappers::VectorBase::equ ( const PetscScalar  a,
const VectorBase V 
)

Assignment *this = a*V.

Definition at line 995 of file petsc_vector_base.cc.

◆ equ() [2/2]

void PETScWrappers::VectorBase::equ ( const PetscScalar  a,
const VectorBase V,
const PetscScalar  b,
const VectorBase W 
)

Assignment *this = a*V + b*W.

This function is deprecated.

Definition at line 1016 of file petsc_vector_base.cc.

◆ ratio()

void PETScWrappers::VectorBase::ratio ( const VectorBase a,
const VectorBase b 
)

Compute the elementwise ratio of the two given vectors, that is let this[i] = a[i]/b[i]. This is useful for example if you want to compute the cellwise ratio of true to estimated error.

This vector is appropriately scaled to hold the result.

If any of the b[i] is zero, the result is undefined. No attempt is made to catch such situations.

Definition at line 1040 of file petsc_vector_base.cc.

◆ write_ascii()

void PETScWrappers::VectorBase::write_ascii ( const PetscViewerFormat  format = PETSC_VIEWER_DEFAULT)

Prints the PETSc vector object values using PETSc internal vector viewer function VecView. The default format prints the vector's contents, including indices of vector elements. For other valid view formats, consult http://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Vec/VecView.html

Definition at line 1051 of file petsc_vector_base.cc.

◆ print()

void PETScWrappers::VectorBase::print ( std::ostream &  out,
const unsigned int  precision = 3,
const bool  scientific = true,
const bool  across = true 
) const

Print to a stream. precision denotes the desired precision with which values shall be printed, scientific whether scientific notation shall be used. If across is true then the vector is printed in a line, while if false then the elements are printed on a separate line each.

Definition at line 1068 of file petsc_vector_base.cc.

◆ swap()

void PETScWrappers::VectorBase::swap ( VectorBase v)

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 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 1115 of file petsc_vector_base.cc.

◆ operator const Vec &()

PETScWrappers::VectorBase::operator const Vec & ( ) const

Conversion operator to gain access to the underlying PETSc type. If you do this, you cut this class off some information it may need, so this conversion operator should only be used if you know what you do. In particular, it should only be used for read-only operations into the vector.

Definition at line 1123 of file petsc_vector_base.cc.

◆ memory_consumption()

std::size_t PETScWrappers::VectorBase::memory_consumption ( ) const

Estimate for the memory consumption (not implemented for this class).

Definition at line 1130 of file petsc_vector_base.cc.

◆ get_mpi_communicator()

virtual const MPI_Comm& PETScWrappers::VectorBase::get_mpi_communicator ( ) const
virtual

Return a reference to the MPI communicator object in use with this object.

Reimplemented in PETScWrappers::MPI::Vector.

◆ do_set_add_operation()

void PETScWrappers::VectorBase::do_set_add_operation ( const size_type  n_elements,
const size_type *  indices,
const PetscScalar *  values,
const bool  add_values 
)
protected

Collective set or add operation: This function is invoked by the collective set and add with the add_values flag set to the corresponding value.

Definition at line 1152 of file petsc_vector_base.cc.

◆ operator=() [2/2]

VectorBase& PETScWrappers::VectorBase::operator= ( const VectorBase )
private

Assignment operator. This is currently not implemented, so it is deliberately left as private (and undefined) to prevent accidental usage.

Friends And Related Function Documentation

◆ internal::VectorReference

friend class internal::VectorReference
friend

Make the reference class a friend.

Definition at line 785 of file petsc_vector_base.h.

◆ swap()

void swap ( VectorBase u,
VectorBase 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 826 of file petsc_vector_base.h.

Member Data Documentation

◆ vector

Vec PETScWrappers::VectorBase::vector
protected

A generic vector object in PETSc. The actual type, a sequential vector, is set in the constructor.

Definition at line 759 of file petsc_vector_base.h.

◆ ghosted

bool PETScWrappers::VectorBase::ghosted
protected

Denotes if this vector has ghost indices associated with it. This means that at least one of the processes in a parallel program has at least one ghost index.

Definition at line 766 of file petsc_vector_base.h.

◆ ghost_indices

IndexSet PETScWrappers::VectorBase::ghost_indices
protected

This vector contains the global indices of the ghost values. The location in this vector denotes the local numbering, which is used in PETSc.

Definition at line 773 of file petsc_vector_base.h.

◆ last_action

VectorOperation::values PETScWrappers::VectorBase::last_action
mutableprotected

Store whether the last action was a write or add operation. This variable is mutable so that the accessor classes can write to it, even though the vector object they refer to is constant.

Definition at line 780 of file petsc_vector_base.h.

◆ attained_ownership

bool PETScWrappers::VectorBase::attained_ownership
protected

Specifies if the vector is the owner of the PETSc Vec. This is true if it got created by this class and determines if it gets destructed in the destructor.

Definition at line 792 of file petsc_vector_base.h.


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