Reference documentation for deal.II version 9.3.3
|
#include <deal.II/lac/petsc_vector.h>
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 MPI_Comm &communicator, const VectorBase &v, const size_type local_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 |
Vector & | operator= (const Vector &v) |
Vector & | operator= (const PetscScalar s) |
template<typename number > | |
Vector & | operator= (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) |
const MPI_Comm & | get_mpi_communicator () const override |
void | print (std::ostream &out, const unsigned int precision=3, const bool scientific=true, const bool across=true) const |
bool | all_zero () const |
void | compress (const VectorOperation::values operation) |
bool | operator== (const VectorBase &v) const |
bool | operator!= (const VectorBase &v) const |
size_type | size () const |
size_type | local_size () const |
size_type | locally_owned_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 |
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 |
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) |
real_type | min () const |
real_type | max () const |
bool | is_non_negative () const |
VectorBase & | operator*= (const PetscScalar factor) |
VectorBase & | operator/= (const PetscScalar factor) |
VectorBase & | operator+= (const VectorBase &V) |
VectorBase & | operator-= (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) |
void | swap (VectorBase &v) |
operator const Vec & () const | |
std::size_t | memory_consumption () const |
Protected Member Functions | |
virtual void | create_vector (const size_type n, const size_type locally_owned_size) |
virtual void | create_vector (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) |
Protected Attributes | |
Vec | vector |
bool | ghosted |
IndexSet | ghost_indices |
VectorOperation::values | last_action |
bool | obtained_ownership |
Private Attributes | |
MPI_Comm | communicator |
Related Functions | |
(Note that these are not member functions.) | |
void | swap (Vector &u, Vector &v) |
void | swap (VectorBase &u, VectorBase &v) |
Subscriptor functionality | |
Classes derived from Subscriptor provide a facility to subscribe to this object. This is mostly used by the SmartPointer class. | |
std::atomic< unsigned int > | counter |
std::map< std::string, unsigned int > | counter_map |
std::vector< std::atomic< bool > * > | validity_pointers |
const std::type_info * | object_info |
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 |
template<class Archive > | |
void | serialize (Archive &ar, const unsigned int version) |
using | map_value_type = decltype(counter_map)::value_type |
using | map_iterator = decltype(counter_map)::iterator |
static std::mutex | mutex |
static ::ExceptionBase & | ExcInUse (int arg1, std::string arg2, std::string arg3) |
static ::ExceptionBase & | ExcNoSubscriber (std::string arg1, std::string arg2) |
void | check_no_subscribers () const noexcept |
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.
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.
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:
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.
Definition at line 156 of file petsc_vector.h.
Declare type for container size.
Definition at line 162 of file petsc_vector.h.
|
inherited |
Declare some of the standard types used in all containers. These types parallel those in the C++
standard libraries vector<...>
class.
Definition at line 258 of file petsc_vector_base.h.
|
inherited |
Definition at line 259 of file petsc_vector_base.h.
|
inherited |
Definition at line 261 of file petsc_vector_base.h.
|
inherited |
Definition at line 262 of file petsc_vector_base.h.
Default constructor. Initialize the vector as empty.
Definition at line 31 of file petsc_parallel_vector.cc.
|
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 42 of file petsc_parallel_vector.cc.
|
explicit |
Copy-constructor from deal.II vectors. Sets the dimension to that of the given vector, and copies all elements.
|
explicit |
Copy-constructor the values from a PETSc wrapper vector class.
Definition at line 52 of file petsc_parallel_vector.cc.
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.
Definition at line 72 of file petsc_parallel_vector.cc.
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 104 of file petsc_parallel_vector.cc.
Copy constructor.
Definition at line 88 of file petsc_parallel_vector.cc.
|
overridevirtual |
Release all memory and return to a state just like after having called the default constructor.
Reimplemented from PETScWrappers::VectorBase.
Definition at line 152 of file petsc_parallel_vector.cc.
Copy the given vector. Resize the present vector if necessary. Also take over the MPI communicator of v
.
Definition at line 115 of file petsc_parallel_vector.cc.
Vector & PETScWrappers::MPI::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.
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.
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 163 of file petsc_parallel_vector.cc.
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 206 of file petsc_parallel_vector.cc.
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.
Definition at line 227 of file petsc_parallel_vector.cc.
Reinit as a vector without ghost elements. See constructor with same signature for more details.
Definition at line 245 of file petsc_parallel_vector.cc.
|
overridevirtual |
Return a reference to the MPI communicator object in use with this vector.
Reimplemented from PETScWrappers::VectorBase.
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.
Definition at line 363 of file petsc_parallel_vector.cc.
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 348 of file petsc_parallel_vector.cc.
|
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 259 of file petsc_parallel_vector.cc.
|
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 277 of file petsc_parallel_vector.cc.
|
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 373 of file petsc_vector_base.cc.
|
inherited |
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 221 of file petsc_vector_base.cc.
|
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 235 of file petsc_vector_base.cc.
|
inherited |
Return the global dimension of the vector.
Definition at line 249 of file petsc_vector_base.cc.
|
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 273 of file petsc_vector_base.cc.
|
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 261 of file petsc_vector_base.cc.
|
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 285 of file petsc_vector_base.cc.
Return whether index
is in the local range or not, see also local_range().
|
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
|
inherited |
Return if the vector contains ghost elements.
|
inherited |
This function only exists for compatibility with the LinearAlgebra::distributed::Vector
class and does nothing: this class implements ghost value updates in a different way that is a better fit with the underlying PETSc vector object.
Provide access to a given element, both read and write.
|
inherited |
Provide read-only access to an element.
Provide access to a given element, both read and write.
Exactly the same as operator().
|
inherited |
Provide read-only access to an element.
Exactly the same as operator().
|
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 298 of file petsc_vector_base.cc.
|
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
indices
and values
arrays must be identical.
|
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
values_begin
as there are iterators between indices_begin
and indices_end
.
|
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 309 of file petsc_vector_base.cc.
|
inherited |
This is a second collective add operation. As a difference, this function takes a deal.II vector of values.
Definition at line 320 of file petsc_vector_base.cc.
|
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 331 of file petsc_vector_base.cc.
|
inherited |
Addition of s
to all components. Note that s
is a scalar and not a vector.
Definition at line 742 of file petsc_vector_base.cc.
|
inherited |
Simple addition of a multiple of a vector, i.e. *this += a*V
.
Definition at line 754 of file petsc_vector_base.cc.
|
inherited |
Multiple addition of scaled vectors, i.e. *this += a*V+b*W
.
Definition at line 766 of file petsc_vector_base.cc.
|
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 340 of file petsc_vector_base.cc.
|
inherited |
Return the square of the \(l_2\)-norm.
Definition at line 433 of file petsc_vector_base.cc.
|
inherited |
Return the mean value of the elements of this vector.
Definition at line 442 of file petsc_vector_base.cc.
|
inherited |
\(l_1\)-norm of the vector. The sum of the absolute values.
Definition at line 493 of file petsc_vector_base.cc.
|
inherited |
\(l_2\)-norm of the vector. The square root of the sum of the squares of the elements.
Definition at line 506 of file petsc_vector_base.cc.
|
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 519 of file petsc_vector_base.cc.
|
inherited |
\(l_\infty\)-norm of the vector. Return the value of the vector element with the maximum absolute value.
Definition at line 561 of file petsc_vector_base.cc.
|
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
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 362 of file petsc_vector_base.cc.
|
inherited |
Return the value of the vector element with the largest negative value.
Definition at line 574 of file petsc_vector_base.cc.
|
inherited |
Return the value of the vector element with the largest positive value.
Definition at line 587 of file petsc_vector_base.cc.
|
inherited |
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 655 of file petsc_vector_base.cc.
|
inherited |
Multiply the entire vector by a fixed factor.
Definition at line 687 of file petsc_vector_base.cc.
|
inherited |
Divide the entire vector by a fixed factor.
Definition at line 701 of file petsc_vector_base.cc.
|
inherited |
Add the given vector to the present one.
Definition at line 718 of file petsc_vector_base.cc.
|
inherited |
Subtract the given vector from the present one.
Definition at line 730 of file petsc_vector_base.cc.
|
inherited |
Scaling and simple vector addition, i.e. this = s(*this)+V
.
Definition at line 785 of file petsc_vector_base.cc.
|
inherited |
Scaling and simple addition, i.e. this = s(*this)+a*V
.
Definition at line 797 of file petsc_vector_base.cc.
|
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 815 of file petsc_vector_base.cc.
|
inherited |
Assignment *this = a*V
.
Definition at line 825 of file petsc_vector_base.cc.
|
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 844 of file petsc_vector_base.cc.
|
inherited |
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 908 of file petsc_vector_base.cc.
|
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 916 of file petsc_vector_base.cc.
|
inherited |
Estimate for the memory consumption (not implemented for this class).
Definition at line 923 of file petsc_vector_base.cc.
|
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 945 of file petsc_vector_base.cc.
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 434 of file petsc_vector.h.
|
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 853 of file petsc_vector_base.h.
|
private |
Copy of the communicator object to be used for this parallel vector.
Definition at line 419 of file petsc_vector.h.
|
protectedinherited |
A generic vector object in PETSc. The actual type, a sequential vector, is set in the constructor.
Definition at line 796 of file petsc_vector_base.h.
|
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 803 of file petsc_vector_base.h.
|
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 810 of file petsc_vector_base.h.
|
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 817 of file petsc_vector_base.h.
|
protectedinherited |
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 destroyed in the destructor.
Definition at line 827 of file petsc_vector_base.h.