Reference documentation for deal.II version GIT relicensing-136-gb80d0be4af 2024-03-18 08:20:02+00:00
\(\newcommand{\dealvcentcolon}{\mathrel{\mathop{:}}}\) \(\newcommand{\dealcoloneq}{\dealvcentcolon\mathrel{\mkern-1.2mu}=}\) \(\newcommand{\jump}[1]{\left[\!\left[ #1 \right]\!\right]}\) \(\newcommand{\average}[1]{\left\{\!\left\{ #1 \right\}\!\right\}}\)
Loading...
Searching...
No Matches
Public Types | Public Member Functions | Static Public Member Functions | Protected Member Functions | Protected Attributes | Private Types | Private Member Functions | Private Attributes | Static Private Attributes | Related Symbols | List of all members
PETScWrappers::MPI::Vector Class Reference

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

Inheritance diagram for PETScWrappers::MPI::Vector:
Inheritance graph
[legend]

Public Types

using size_type = types::global_dof_index
 
using value_type = PetscScalar
 
using real_type = PetscReal
 
using reference = internal::VectorReference
 
using const_reference = const internal::VectorReference
 

Public Member Functions

 Vector ()
 
 Vector (const MPI_Comm communicator, const size_type n, const size_type locally_owned_size)
 
template<typename Number >
 Vector (const MPI_Comm communicator, const ::Vector< Number > &v, const size_type locally_owned_size)
 
 Vector (const IndexSet &local, const IndexSet &ghost, const MPI_Comm communicator)
 
 Vector (const IndexSet &local, const MPI_Comm communicator)
 
 Vector (const Vector &v)
 
virtual void clear () override
 
Vectoroperator= (const Vector &v)
 
Vectoroperator= (const PetscScalar s)
 
template<typename number >
Vectoroperator= (const ::Vector< number > &v)
 
void reinit (const MPI_Comm communicator, const size_type N, const size_type locally_owned_size, const bool omit_zeroing_entries=false)
 
void reinit (const Vector &v, const bool omit_zeroing_entries=false)
 
void reinit (const IndexSet &local, const IndexSet &ghost, const MPI_Comm communicator)
 
void reinit (const IndexSet &local, const MPI_Comm communicator)
 
void reinit (const std::shared_ptr< const Utilities::MPI::Partitioner > &partitioner, const bool make_ghosted=true)
 
void print (std::ostream &out, const unsigned int precision=3, const bool scientific=true, const bool across=true) const
 
bool all_zero () const
 
 VectorBase ()
 
 VectorBase (const VectorBase &v)
 
 VectorBase (const Vec &v)
 
void reinit (Vec v)
 
void compress (const VectorOperation::values operation)
 
bool operator== (const VectorBase &v) const
 
bool operator!= (const VectorBase &v) const
 
size_type size () const override
 
size_type locally_owned_size () const
 
std::pair< size_type, size_typelocal_range () const
 
bool in_local_range (const size_type index) const
 
IndexSet locally_owned_elements () const
 
bool has_ghost_elements () const
 
const IndexSetghost_elements () const
 
void update_ghost_values () 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
 
virtual void extract_subvector_to (const ArrayView< const types::global_dof_index > &indices, ArrayView< PetscScalar > &elements) const override
 
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)
 
void add (const PetscScalar s)
 
void add (const PetscScalar a, const VectorBase &V)
 
void add (const PetscScalar a, const VectorBase &V, const PetscScalar b, const VectorBase &W)
 
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)
 
VectorBaseoperator*= (const PetscScalar factor)
 
VectorBaseoperator/= (const PetscScalar factor)
 
VectorBaseoperator+= (const VectorBase &V)
 
VectorBaseoperator-= (const VectorBase &V)
 
void sadd (const PetscScalar s, const VectorBase &V)
 
void sadd (const PetscScalar s, const PetscScalar a, const VectorBase &V)
 
void scale (const VectorBase &scaling_factors)
 
void equ (const PetscScalar a, const VectorBase &V)
 
void write_ascii (const PetscViewerFormat format=PETSC_VIEWER_DEFAULT)
 
template<class Archive >
void save (Archive &ar, const unsigned int version) const
 
template<class Archive >
void load (Archive &ar, const unsigned int version)
 
template<class Archive >
void serialize (Archive &archive, const unsigned int version)
 
void swap (VectorBase &v) noexcept
 
 operator const Vec & () const
 
Vec & petsc_vector ()
 
std::size_t memory_consumption () const
 
MPI_Comm get_mpi_communicator () const
 
Subscriptor functionality

Classes derived from Subscriptor provide a facility to subscribe to this object. This is mostly used by the SmartPointer class.

void subscribe (std::atomic< bool > *const validity, const std::string &identifier="") const
 
void unsubscribe (std::atomic< bool > *const validity, const std::string &identifier="") const
 
unsigned int n_subscriptions () const
 
template<typename StreamType >
void list_subscribers (StreamType &stream) const
 
void list_subscribers () const
 

Static Public Member Functions

static ::ExceptionBaseExcInUse (int arg1, std::string arg2, std::string arg3)
 
static ::ExceptionBaseExcNoSubscriber (std::string arg1, std::string arg2)
 

Protected Member Functions

virtual void create_vector (const MPI_Comm comm, const size_type n, const size_type locally_owned_size)
 
virtual void create_vector (const MPI_Comm comm, const size_type n, const size_type locally_owned_size, const IndexSet &ghostnodes)
 
void do_set_add_operation (const size_type n_elements, const size_type *indices, const PetscScalar *values, const bool add_values)
 
void determine_ghost_indices ()
 

Protected Attributes

Vec vector
 
bool ghosted
 
IndexSet ghost_indices
 
VectorOperation::values last_action
 

Private Types

using map_value_type = decltype(counter_map)::value_type
 
using map_iterator = decltype(counter_map)::iterator
 

Private Member Functions

void check_no_subscribers () const noexcept
 

Private Attributes

std::atomic< unsigned intcounter
 
std::map< std::string, unsigned intcounter_map
 
std::vector< std::atomic< bool > * > validity_pointers
 
const std::type_info * object_info
 

Static Private Attributes

static std::mutex mutex
 

Related Symbols

(Note that these are not member symbols.)

void swap (Vector &u, Vector &v)
 
void swap (VectorBase &u, VectorBase &v) noexcept
 

Detailed Description

Implementation of a parallel vector class based on PETSC and using MPI communication to synchronize distributed operations. All the functionality is actually in the base class, except for the calls to generate a parallel 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.

Parallel communication model

The parallel functionality of PETSc is built on top of the Message Passing Interface (MPI). MPI's communication model is built on collective communications: if one process wants something from another, that other process has to be willing to accept this communication. A process cannot query data from another process by calling a remote function, without that other process expecting such a transaction. The consequence is that most of the operations in the base class of this class have to be called collectively. For example, if you want to compute the l2 norm of a parallel vector, all processes across which this vector is shared have to call the l2_norm function. If you don't do this, but instead only call the l2_norm function on one process, then the following happens: This one process will call one of the collective MPI functions and wait for all the other processes to join in on this. Since the other processes don't call this function, you will either get a time-out on the first process, or, worse, by the time the next a call to a PETSc function generates an MPI message on the other processes, you will get a cryptic message that only a subset of processes attempted a communication. These bugs can be very hard to figure out, unless you are well-acquainted with the communication model of MPI, and know which functions may generate MPI messages.

One particular case, where an MPI message may be generated unexpectedly is discussed below.

Accessing individual elements of a vector

PETSc does allow read access to individual elements of a vector, but in the distributed case only to elements that are stored locally. We implement this through calls like d=vec(i). However, if you access an element outside the locally stored range, an exception is generated.

In contrast to read access, PETSc (and the respective deal.II wrapper classes) allow to write (or add) to individual elements of vectors, even if they are stored on a different process. You can do this writing, for example, vec(i)=d or vec(i)+=d, or similar operations. There is one catch, however, that may lead to very confusing error messages: PETSc requires application programs to call the compress() function when they switch from adding, to elements to writing to elements. The reasoning is that all processes might accumulate addition operations to elements, even if multiple processes write to the same elements. By the time we call compress() the next time, all these additions are executed. However, if one process adds to an element, and another overwrites to it, the order of execution would yield non-deterministic behavior if we don't make sure that a synchronization with compress() happens in between.

In order to make sure these calls to compress() happen at the appropriate time, the deal.II wrappers keep a state variable that store which is the presently allowed operation: additions or writes. If it encounters an operation of the opposite kind, it calls compress() and flips the state. This can sometimes lead to very confusing behavior, in code that may for example look like this:

...
// do some write operations on the vector
for (unsigned int i=0; i<vector.size(); ++i)
vector(i) = i;
// do some additions to vector elements, but only for some elements
for (unsigned int i=0; i<vector.size(); ++i)
if (some_condition(i) == true)
vector(i) += 1;
// do another collective operation
const double norm = vector.l2_norm();

This code can run into trouble: by the time we see the first addition operation, we need to flush the overwrite buffers for the vector, and the deal.II library will do so by calling compress(). However, it will only do so for all processes that actually do an addition – if the condition is never true for one of the processes, then this one will not get to the actual compress() call, whereas all the other ones do. This gets us into trouble, since all the other processes hang in the call to flush the write buffers, while the one other process advances to the call to compute the l2 norm. At this time, you will get an error that some operation was attempted by only a subset of processes. This behavior may seem surprising, unless you know that write/addition operations on single elements may trigger this behavior.

The problem described here may be avoided by placing additional calls to compress(), or making sure that all processes do the same type of operations at the same time, for example by placing zero additions if necessary.

See also
vectors with ghost elements

Definition at line 157 of file petsc_vector.h.

Member Typedef Documentation

◆ size_type

Declare type for container size.

Definition at line 163 of file petsc_vector.h.

◆ 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 259 of file petsc_vector_base.h.

◆ real_type

Definition at line 260 of file petsc_vector_base.h.

◆ reference

using PETScWrappers::VectorBase::reference = internal::VectorReference
inherited

Definition at line 262 of file petsc_vector_base.h.

◆ const_reference

using PETScWrappers::VectorBase::const_reference = const internal::VectorReference
inherited

Definition at line 263 of file petsc_vector_base.h.

◆ map_value_type

The data type used in counter_map.

Definition at line 229 of file subscriptor.h.

◆ map_iterator

using Subscriptor::map_iterator = decltype(counter_map)::iterator
privateinherited

The iterator type used in counter_map.

Definition at line 234 of file subscriptor.h.

Constructor & Destructor Documentation

◆ Vector() [1/6]

Vector< Number >::Vector ( )

Default constructor. Initialize the vector as empty.

Definition at line 38 of file petsc_parallel_vector.cc.

◆ Vector() [2/6]

Vector< Number >::Vector ( const MPI_Comm  communicator,
const size_type  n,
const size_type  locally_owned_size 
)
explicit

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

  • locally_owned_size denotes the size of the chunk that shall be stored on the present process.
  • communicator denotes the MPI communicator over which the different parts of the vector shall communicate

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 48 of file petsc_parallel_vector.cc.

◆ Vector() [3/6]

template<typename Number >
PETScWrappers::MPI::Vector::Vector ( const MPI_Comm  communicator,
const ::Vector< Number > &  v,
const size_type  locally_owned_size 
)
explicit

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

  • locally_owned_size denotes the size of the chunk that shall be stored on the present process.
  • communicator denotes the MPI communicator over which the different parts of the vector shall communicate

◆ Vector() [4/6]

Vector< Number >::Vector ( const IndexSet local,
const IndexSet ghost,
const MPI_Comm  communicator 
)

Construct a new parallel ghosted PETSc vector from IndexSets.

Note that local must be ascending and 1:1, see IndexSet::is_ascending_and_one_to_one(). In particular, the DoFs in local need to be contiguous, meaning you can only create vectors from a DoFHandler with several finite element components if they are not reordered by component (use a PETScWrappers::BlockVector otherwise). The global size of the vector is determined by local.size(). The global indices in ghost are supplied as ghost indices so that they can be read locally.

Note that the ghost IndexSet may be empty and that any indices already contained in local are ignored during construction. That way, the ghost parameter can equal the set of locally relevant degrees of freedom, see step-32.

Note
This operation always creates a ghosted vector, which is considered read-only.
See also
vectors with ghost elements

Definition at line 57 of file petsc_parallel_vector.cc.

◆ Vector() [5/6]

Vector< Number >::Vector ( const IndexSet local,
const MPI_Comm  communicator 
)
explicit

Construct a new parallel PETSc vector without ghost elements from an IndexSet.

Note that local must be ascending and 1:1, see IndexSet::is_ascending_and_one_to_one(). In particular, the DoFs in local need to be contiguous, meaning you can only create vectors from a DoFHandler with several finite element components if they are not reordered by component (use a PETScWrappers::BlockVector otherwise).

Definition at line 93 of file petsc_parallel_vector.cc.

◆ Vector() [6/6]

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

Copy constructor.

Definition at line 75 of file petsc_parallel_vector.cc.

Member Function Documentation

◆ clear()

void Vector< Number >::clear ( )
overridevirtual

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

Reimplemented from PETScWrappers::VectorBase.

Definition at line 145 of file petsc_parallel_vector.cc.

◆ operator=() [1/3]

Vector & Vector< Number >::operator= ( const Vector v)

Copy the given vector. Resize the present vector if necessary. Also take over the MPI communicator of v.

Definition at line 103 of file petsc_parallel_vector.cc.

◆ operator=() [2/3]

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.

◆ operator=() [3/3]

template<typename number >
Vector & PETScWrappers::MPI::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.

Contrary to the case of sequential vectors, this operators requires that the present vector already has the correct size, since we need to have a partition and a communicator present which we otherwise can't get from the source vector.

◆ reinit() [1/6]

void Vector< Number >::reinit ( const MPI_Comm  communicator,
const size_type  N,
const size_type  locally_owned_size,
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

locally_owned_size denotes how many of the N values shall be stored locally on the present process. for less data.

communicator denotes the MPI communicator henceforth to be used for this vector.

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

Definition at line 155 of file petsc_parallel_vector.cc.

◆ reinit() [2/6]

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

Change the dimension to that of the vector v, and also take over the partitioning into local sizes as well as the MPI communicator. 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(), v.locally_owned_size(), omit_zeroing_entries).

Definition at line 196 of file petsc_parallel_vector.cc.

◆ reinit() [3/6]

void Vector< Number >::reinit ( const IndexSet local,
const IndexSet ghost,
const MPI_Comm  communicator 
)

Reinit as a vector with ghost elements. See the constructor with same signature for more details.

See also
vectors with ghost elements

Definition at line 219 of file petsc_parallel_vector.cc.

◆ reinit() [4/6]

void Vector< Number >::reinit ( const IndexSet local,
const MPI_Comm  communicator 
)

Reinit as a vector without ghost elements. See constructor with same signature for more details.

See also
vectors with ghost elements

Definition at line 235 of file petsc_parallel_vector.cc.

◆ reinit() [5/6]

void Vector< Number >::reinit ( const std::shared_ptr< const Utilities::MPI::Partitioner > &  partitioner,
const bool  make_ghosted = true 
)

Initialize the vector given to the parallel partitioning described in partitioner.

You can decide whether your vector will contain ghost elements with make_ghosted.

Definition at line 246 of file petsc_parallel_vector.cc.

◆ print()

void Vector< Number >::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.

Note
This function overloads the one in the base class to ensure that the right thing happens for parallel vectors that are distributed across processors.

Definition at line 360 of file petsc_parallel_vector.cc.

◆ all_zero()

bool Vector< Number >::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.

Note
This function overloads the one in the base class to make this a collective operation.

Definition at line 348 of file petsc_parallel_vector.cc.

◆ create_vector() [1/2]

void Vector< Number >::create_vector ( const MPI_Comm  comm,
const size_type  n,
const size_type  locally_owned_size 
)
protectedvirtual

Create a vector of length n. For this class, we create a parallel vector. n denotes the total size of the vector to be created. locally_owned_size denotes how many of these elements shall be stored locally.

Definition at line 269 of file petsc_parallel_vector.cc.

◆ create_vector() [2/2]

void Vector< Number >::create_vector ( const MPI_Comm  comm,
const size_type  n,
const size_type  locally_owned_size,
const IndexSet ghostnodes 
)
protectedvirtual

Create a vector of global length n, local size locally_owned_size and with the specified ghost indices. Note that you need to call update_ghost_values() before accessing those.

Definition at line 289 of file petsc_parallel_vector.cc.

◆ VectorBase() [1/3]

PETScWrappers::VectorBase::VectorBase ( )

Import VectorBase constructors, including from a PETSc Vec object.

Definition at line 269 of file petsc_vector_base.cc.

◆ VectorBase() [2/3]

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

Import VectorBase constructors, including from a PETSc Vec object.

Definition at line 275 of file petsc_vector_base.cc.

◆ VectorBase() [3/3]

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

Import VectorBase constructors, including from a PETSc Vec object.

Definition at line 281 of file petsc_vector_base.cc.

◆ reinit() [6/6]

void PETScWrappers::VectorBase::reinit ( Vec  v)

This method associates the PETSc Vec to the instance of the class. This is particularly useful when performing PETSc to Deal.II operations since it allows to reuse the Deal.II VectorBase and the PETSc Vec without incurring in memory copies.

Definition at line 337 of file petsc_vector_base.cc.

◆ compress()

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

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

◆ operator==()

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

◆ operator!=()

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

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

◆ size()

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

Return the global dimension of the vector.

Implements ReadVector< PetscScalar >.

Definition at line 427 of file petsc_vector_base.cc.

◆ locally_owned_size()

VectorBase::size_type PETScWrappers::VectorBase::locally_owned_size ( ) const
inherited

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() or locally_owned_elements().

Definition at line 439 of file petsc_vector_base.cc.

◆ local_range()

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

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=locally_owned_size().

Definition at line 451 of file petsc_vector_base.cc.

◆ in_local_range()

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

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

◆ locally_owned_elements()

IndexSet PETScWrappers::VectorBase::locally_owned_elements ( ) const
inherited

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())
IndexSet complete_index_set(const IndexSet::size_type N)
Definition index_set.h:1193

◆ has_ghost_elements()

bool PETScWrappers::VectorBase::has_ghost_elements ( ) const
inherited

Return if the vector contains ghost elements.

See also
vectors with ghost elements

◆ ghost_elements()

const IndexSet & PETScWrappers::VectorBase::ghost_elements ( ) const
inherited

Return the IndexSet of ghost elements.

◆ update_ghost_values()

void PETScWrappers::VectorBase::update_ghost_values ( ) const
inherited

Update ghosted elements.

◆ operator()() [1/2]

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

◆ operator()() [2/2]

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

Provide read-only access to an element.

◆ operator[]() [1/2]

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
inherited

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 
)
inherited

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

◆ extract_subvector_to() [1/3]

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

Instead of getting individual elements of a vector via operator(), this function allows getting 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.

If the current vector is called v, then this function is the equivalent to the code

for (unsigned int i=0; i<indices.size(); ++i)
values[i] = v[indices[i]];
Precondition
The sizes of the indices and values arrays must be identical.

◆ extract_subvector_to() [2/3]

virtual void PETScWrappers::VectorBase::extract_subvector_to ( const ArrayView< const types::global_dof_index > &  indices,
ArrayView< PetscScalar > &  elements 
) const
overridevirtualinherited

Extract a range of elements all at once.

Implements ReadVector< PetscScalar >.

◆ extract_subvector_to() [3/3]

void PETScWrappers::VectorBase::extract_subvector_to ( const ForwardIterator  indices_begin,
const ForwardIterator  indices_end,
OutputIterator  values_begin 
) const
inherited

Instead of getting individual elements of a vector via operator(), this function allows getting a whole set of elements at once. In contrast to the previous function, this function obtains the indices of the elements by dereferencing all elements of the iterator range provided by the first two arguments, and puts the vector values into memory locations obtained by dereferencing a range of iterators starting at the location pointed to by the third argument.

If the current vector is called v, then this function is the equivalent to the code

Precondition
It must be possible to write into as many memory locations starting at values_begin as there are iterators between indices_begin and indices_end.

◆ add() [1/6]

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

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

◆ add() [2/6]

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

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

Definition at line 486 of file petsc_vector_base.cc.

◆ add() [3/6]

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

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

◆ add() [4/6]

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

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

Definition at line 855 of file petsc_vector_base.cc.

◆ add() [5/6]

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

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

Definition at line 867 of file petsc_vector_base.cc.

◆ add() [6/6]

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

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

Definition at line 879 of file petsc_vector_base.cc.

◆ operator*()

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

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

◆ norm_sqr()

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

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

Definition at line 604 of file petsc_vector_base.cc.

◆ mean_value()

PetscScalar PETScWrappers::VectorBase::mean_value ( ) const
inherited

Return the mean value of the elements of this vector.

Definition at line 613 of file petsc_vector_base.cc.

◆ l1_norm()

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

\(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 664 of file petsc_vector_base.cc.

◆ l2_norm()

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

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

Definition at line 677 of file petsc_vector_base.cc.

◆ lp_norm()

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

\(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 690 of file petsc_vector_base.cc.

◆ linfty_norm()

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

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

Definition at line 732 of file petsc_vector_base.cc.

◆ add_and_dot()

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

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;
void add(const std::vector< size_type > &indices, const std::vector< PetscScalar > &values)

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.

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

◆ operator*=()

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

Multiply the entire vector by a fixed factor.

Definition at line 800 of file petsc_vector_base.cc.

◆ operator/=()

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

Divide the entire vector by a fixed factor.

Definition at line 814 of file petsc_vector_base.cc.

◆ operator+=()

Add the given vector to the present one.

Definition at line 831 of file petsc_vector_base.cc.

◆ operator-=()

Subtract the given vector from the present one.

Definition at line 843 of file petsc_vector_base.cc.

◆ sadd() [1/2]

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

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

Definition at line 898 of file petsc_vector_base.cc.

◆ sadd() [2/2]

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

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

Definition at line 910 of file petsc_vector_base.cc.

◆ scale()

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

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

◆ equ()

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

Assignment *this = a*V.

Definition at line 938 of file petsc_vector_base.cc.

◆ write_ascii()

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

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

◆ save()

template<class Archive >
void PETScWrappers::VectorBase::save ( Archive ar,
const unsigned int  version 
) const
inherited

Write the data of this object to a stream for the purpose of serialization using the BOOST serialization library.

Note
Each processor only serializes its own locally owned values.

◆ load()

template<class Archive >
void PETScWrappers::VectorBase::load ( Archive ar,
const unsigned int  version 
)
inherited

Read the data of this object from a stream for the purpose of serialization using the BOOST serialization library.

◆ serialize()

template<class Archive >
void PETScWrappers::VectorBase::serialize ( Archive archive,
const unsigned int  version 
)
inherited

Write and read the data of this object from a stream for the purpose of serialization using the BOOST serialization library.

◆ swap()

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

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

◆ operator const Vec &()

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

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

◆ petsc_vector()

Vec & PETScWrappers::VectorBase::petsc_vector ( )
inherited

Return a reference to the underlying PETSc type. It can be used to modify the underlying data, so use it only when you know what you are doing.

Definition at line 1036 of file petsc_vector_base.cc.

◆ memory_consumption()

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

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

Definition at line 1043 of file petsc_vector_base.cc.

◆ get_mpi_communicator()

MPI_Comm PETScWrappers::VectorBase::get_mpi_communicator ( ) const
inherited

Return the underlying MPI communicator.

◆ 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 
)
protectedinherited

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

◆ determine_ghost_indices()

void PETScWrappers::VectorBase::determine_ghost_indices ( )
protectedinherited

Determine ghost indices from the internal PETSc Vec

Definition at line 242 of file petsc_vector_base.cc.

◆ subscribe()

void Subscriptor::subscribe ( std::atomic< bool > *const  validity,
const std::string &  identifier = "" 
) const
inherited

Subscribes a user of the object by storing the pointer validity. The subscriber may be identified by text supplied as identifier.

Definition at line 135 of file subscriptor.cc.

◆ unsubscribe()

void Subscriptor::unsubscribe ( std::atomic< bool > *const  validity,
const std::string &  identifier = "" 
) const
inherited

Unsubscribes a user from the object.

Note
The identifier and the validity pointer must be the same as the one supplied to subscribe().

Definition at line 155 of file subscriptor.cc.

◆ n_subscriptions()

unsigned int Subscriptor::n_subscriptions ( ) const
inlineinherited

Return the present number of subscriptions to this object. This allows to use this class for reference counted lifetime determination where the last one to unsubscribe also deletes the object.

Definition at line 300 of file subscriptor.h.

◆ list_subscribers() [1/2]

template<typename StreamType >
void Subscriptor::list_subscribers ( StreamType &  stream) const
inlineinherited

List the subscribers to the input stream.

Definition at line 317 of file subscriptor.h.

◆ list_subscribers() [2/2]

void Subscriptor::list_subscribers ( ) const
inherited

List the subscribers to deallog.

Definition at line 203 of file subscriptor.cc.

◆ check_no_subscribers()

void Subscriptor::check_no_subscribers ( ) const
privatenoexceptinherited

Check that there are no objects subscribing to this object. If this check passes then it is safe to destroy the current object. It this check fails then this function will either abort or print an error message to deallog (by using the AssertNothrow mechanism), but will not throw an exception.

Note
Since this function is just a consistency check it does nothing in release mode.
If this function is called when there is an uncaught exception then, rather than aborting, this function prints an error message to the standard error stream and returns.

Definition at line 52 of file subscriptor.cc.

Friends And Related Symbol Documentation

◆ swap() [1/2]

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.

Definition at line 425 of file petsc_vector.h.

◆ swap() [2/2]

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.

Definition at line 868 of file petsc_vector_base.h.

Member Data Documentation

◆ vector

Vec PETScWrappers::VectorBase::vector
protectedinherited

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

Definition at line 812 of file petsc_vector_base.h.

◆ ghosted

bool PETScWrappers::VectorBase::ghosted
protectedinherited

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 819 of file petsc_vector_base.h.

◆ ghost_indices

IndexSet PETScWrappers::VectorBase::ghost_indices
protectedinherited

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 826 of file petsc_vector_base.h.

◆ last_action

VectorOperation::values PETScWrappers::VectorBase::last_action
mutableprotectedinherited

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 833 of file petsc_vector_base.h.

◆ counter

std::atomic<unsigned int> Subscriptor::counter
mutableprivateinherited

Store the number of objects which subscribed to this object. Initially, this number is zero, and upon destruction it shall be zero again (i.e. all objects which subscribed should have unsubscribed again).

The creator (and owner) of an object is counted in the map below if HE manages to supply identification.

We use the mutable keyword in order to allow subscription to constant objects also.

This counter may be read from and written to concurrently in multithreaded code: hence we use the std::atomic class template.

Definition at line 218 of file subscriptor.h.

◆ counter_map

std::map<std::string, unsigned int> Subscriptor::counter_map
mutableprivateinherited

In this map, we count subscriptions for each different identification string supplied to subscribe().

Definition at line 224 of file subscriptor.h.

◆ validity_pointers

std::vector<std::atomic<bool> *> Subscriptor::validity_pointers
mutableprivateinherited

In this vector, we store pointers to the validity bool in the SmartPointer objects that subscribe to this class.

Definition at line 240 of file subscriptor.h.

◆ object_info

const std::type_info* Subscriptor::object_info
mutableprivateinherited

Pointer to the typeinfo object of this object, from which we can later deduce the class name. Since this information on the derived class is neither available in the destructor, nor in the constructor, we obtain it in between and store it here.

Definition at line 248 of file subscriptor.h.

◆ mutex

std::mutex Subscriptor::mutex
staticprivateinherited

A mutex used to ensure data consistency when accessing the mutable members of this class. This lock is used in the subscribe() and unsubscribe() functions, as well as in list_subscribers().

Definition at line 271 of file subscriptor.h.


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