Reference documentation for deal.II version 9.6.0
|
#include <deal.II/lac/petsc_vector_base.h>
Public Types | |
using | value_type = PetscScalar |
using | real_type = PetscReal |
using | size_type = types::global_dof_index |
using | reference = internal::VectorReference |
using | const_reference = const internal::VectorReference |
Public Member Functions | |
VectorBase () | |
VectorBase (const VectorBase &v) | |
VectorBase (const Vec &v) | |
virtual | ~VectorBase () override |
virtual void | clear () |
void | compress (const VectorOperation::values operation) |
VectorBase & | operator= (const VectorBase &) |
VectorBase & | operator= (const PetscScalar s) |
void | reinit (Vec v) |
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_type > | local_range () const |
bool | in_local_range (const size_type index) const |
IndexSet | locally_owned_elements () const |
bool | has_ghost_elements () const |
const IndexSet & | 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 |
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) |
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) |
bool | all_zero () const |
VectorBase & | operator*= (const PetscScalar factor) |
VectorBase & | operator/= (const PetscScalar factor) |
VectorBase & | operator+= (const VectorBase &V) |
VectorBase & | operator-= (const VectorBase &V) |
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) |
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 | print (std::ostream &out, const unsigned int precision=3, const bool scientific=true, const bool across=true) const |
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 ::ExceptionBase & | ExcInUse (int arg1, std::string arg2, std::string arg3) |
static ::ExceptionBase & | ExcNoSubscriber (std::string arg1, std::string arg2) |
Protected Member Functions | |
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 int > | counter |
std::map< std::string, unsigned int > | counter_map |
std::vector< std::atomic< bool > * > | validity_pointers |
const std::type_info * | object_info |
Static Private Attributes | |
static std::mutex | mutex |
Friends | |
class | internal::VectorReference |
Related Symbols | |
(Note that these are not member symbols.) | |
void | swap (VectorBase &u, VectorBase &v) noexcept |
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. As such, this class can also be used as a deal.II-compatible wrapper for a PETSc Vec
object of any type. 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.
Definition at line 252 of file petsc_vector_base.h.
using PETScWrappers::VectorBase::value_type = PetscScalar |
Declare some of the standard types used in all containers. These types parallel those in the C++
standard libraries vector<...>
class.
Definition at line 260 of file petsc_vector_base.h.
using PETScWrappers::VectorBase::real_type = PetscReal |
Definition at line 261 of file petsc_vector_base.h.
Definition at line 262 of file petsc_vector_base.h.
using PETScWrappers::VectorBase::reference = internal::VectorReference |
Definition at line 263 of file petsc_vector_base.h.
using PETScWrappers::VectorBase::const_reference = const internal::VectorReference |
Definition at line 264 of file petsc_vector_base.h.
|
privateinherited |
The data type used in counter_map.
Definition at line 229 of file subscriptor.h.
|
privateinherited |
The iterator type used in counter_map.
Definition at line 234 of file subscriptor.h.
PETScWrappers::VectorBase::VectorBase | ( | ) |
Default constructor. It doesn't do anything, derived classes will have to initialize the data.
Definition at line 122 of file petsc_vector_base.cc.
PETScWrappers::VectorBase::VectorBase | ( | const VectorBase & | v | ) |
Copy constructor. Sets the dimension to that of the given vector, and copies all elements.
Definition at line 130 of file petsc_vector_base.cc.
|
explicit |
Initialize a Vector from a PETSc Vec object. Note that we do not copy the vector.
Definition at line 145 of file petsc_vector_base.cc.
|
overridevirtual |
Destructor.
Definition at line 160 of file petsc_vector_base.cc.
|
virtual |
Release all memory and return to a state just like after having called the default constructor.
Reimplemented in PETScWrappers::MPI::Vector.
Definition at line 345 of file petsc_vector_base.cc.
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 540 of file petsc_vector_base.cc.
VectorBase & PETScWrappers::VectorBase::operator= | ( | const VectorBase & | v | ) |
The copy assignment operator.
Definition at line 358 of file petsc_vector_base.cc.
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 371 of file petsc_vector_base.cc.
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 170 of file petsc_vector_base.cc.
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 399 of file petsc_vector_base.cc.
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 413 of file petsc_vector_base.cc.
|
overridevirtual |
Return the global dimension of the vector.
Implements ReadVector< PetscScalar >.
Definition at line 427 of file petsc_vector_base.cc.
VectorBase::size_type PETScWrappers::VectorBase::locally_owned_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() or locally_owned_elements().
Definition at line 439 of file petsc_vector_base.cc.
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=locally_owned_size()
.
Definition at line 451 of file petsc_vector_base.cc.
Return whether index
is in the local range or not, see also local_range().
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
bool PETScWrappers::VectorBase::has_ghost_elements | ( | ) | const |
Return if the vector contains ghost elements.
const IndexSet & PETScWrappers::VectorBase::ghost_elements | ( | ) | const |
Return the IndexSet of ghost elements.
void PETScWrappers::VectorBase::update_ghost_values | ( | ) | const |
Update ghosted elements.
Provide access to a given element, both read and write.
PetscScalar PETScWrappers::VectorBase::operator() | ( | const size_type | index | ) | const |
Provide read-only access to an element.
Provide access to a given element, both read and write.
Exactly the same as operator().
PetscScalar PETScWrappers::VectorBase::operator[] | ( | const size_type | index | ) | const |
Provide read-only access to an element.
Exactly the same as operator().
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 464 of file petsc_vector_base.cc.
void PETScWrappers::VectorBase::extract_subvector_to | ( | const std::vector< size_type > & | indices, |
std::vector< PetscScalar > & | values ) const |
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.
|
overridevirtual |
Extract a range of elements all at once.
Implements ReadVector< PetscScalar >.
void PETScWrappers::VectorBase::extract_subvector_to | ( | const ForwardIterator | indices_begin, |
const ForwardIterator | indices_end, | ||
OutputIterator | values_begin ) const |
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
. 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 475 of file petsc_vector_base.cc.
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 486 of file petsc_vector_base.cc.
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 497 of file petsc_vector_base.cc.
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 507 of file petsc_vector_base.cc.
VectorBase::real_type PETScWrappers::VectorBase::norm_sqr | ( | ) | const |
Return the square of the \(l_2\)-norm.
Definition at line 604 of file petsc_vector_base.cc.
PetscScalar PETScWrappers::VectorBase::mean_value | ( | ) | const |
Return the mean value of the elements of this vector.
Definition at line 613 of file petsc_vector_base.cc.
VectorBase::real_type PETScWrappers::VectorBase::l1_norm | ( | ) | const |
\(l_1\)-norm of the vector. The sum of the absolute values.
Definition at line 664 of file petsc_vector_base.cc.
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 677 of file petsc_vector_base.cc.
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 690 of file petsc_vector_base.cc.
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 732 of file petsc_vector_base.cc.
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
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.
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 745 of file petsc_vector_base.cc.
VectorBase & PETScWrappers::VectorBase::operator*= | ( | const PetscScalar | factor | ) |
Multiply the entire vector by a fixed factor.
Definition at line 800 of file petsc_vector_base.cc.
VectorBase & PETScWrappers::VectorBase::operator/= | ( | const PetscScalar | factor | ) |
Divide the entire vector by a fixed factor.
Definition at line 814 of file petsc_vector_base.cc.
VectorBase & PETScWrappers::VectorBase::operator+= | ( | const VectorBase & | V | ) |
Add the given vector to the present one.
Definition at line 831 of file petsc_vector_base.cc.
VectorBase & PETScWrappers::VectorBase::operator-= | ( | const VectorBase & | V | ) |
Subtract the given vector from the present one.
Definition at line 843 of file petsc_vector_base.cc.
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 855 of file petsc_vector_base.cc.
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 867 of file petsc_vector_base.cc.
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 879 of file petsc_vector_base.cc.
void PETScWrappers::VectorBase::sadd | ( | const PetscScalar | s, |
const VectorBase & | V ) |
Scaling and simple vector addition, i.e. this = s(*this)+V
.
Definition at line 898 of file petsc_vector_base.cc.
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 910 of file petsc_vector_base.cc.
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 928 of file petsc_vector_base.cc.
void PETScWrappers::VectorBase::equ | ( | const PetscScalar | a, |
const VectorBase & | V ) |
Assignment *this = a*V
.
Definition at line 938 of file petsc_vector_base.cc.
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 952 of file petsc_vector_base.cc.
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 972 of file petsc_vector_base.cc.
void PETScWrappers::VectorBase::save | ( | Archive & | ar, |
const unsigned int | version ) const |
Write the data of this object to a stream for the purpose of serialization using the BOOST serialization library.
void PETScWrappers::VectorBase::load | ( | Archive & | ar, |
const unsigned int | version ) |
Read the data of this object from a stream for the purpose of serialization using the BOOST serialization library.
void PETScWrappers::VectorBase::serialize | ( | Archive & | archive, |
const unsigned int | version ) |
Write and read the data of this object from a stream for the purpose of serialization using the BOOST serialization library.
|
noexcept |
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 1019 of file petsc_vector_base.cc.
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 1031 of file petsc_vector_base.cc.
Vec & PETScWrappers::VectorBase::petsc_vector | ( | ) |
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 1038 of file petsc_vector_base.cc.
std::size_t PETScWrappers::VectorBase::memory_consumption | ( | ) | const |
Estimate for the memory consumption (not implemented for this class).
Definition at line 1045 of file petsc_vector_base.cc.
MPI_Comm PETScWrappers::VectorBase::get_mpi_communicator | ( | ) | const |
Return the underlying MPI communicator.
|
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 1067 of file petsc_vector_base.cc.
|
protected |
Determine ghost indices from the internal PETSc Vec
Definition at line 242 of file petsc_vector_base.cc.
|
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.
|
inherited |
Unsubscribes a user from the object.
identifier
and the validity
pointer must be the same as the one supplied to subscribe(). Definition at line 155 of file subscriptor.cc.
|
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.
|
inlineinherited |
List the subscribers to the input stream
.
Definition at line 317 of file subscriptor.h.
|
inherited |
List the subscribers to deallog
.
Definition at line 203 of file subscriptor.cc.
|
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.
Definition at line 52 of file subscriptor.cc.
|
friend |
Definition at line 837 of file petsc_vector_base.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 869 of file petsc_vector_base.h.
|
protected |
A generic vector object in PETSc. The actual type, a sequential vector, is set in the constructor.
Definition at line 813 of file petsc_vector_base.h.
|
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 820 of file petsc_vector_base.h.
|
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 827 of file petsc_vector_base.h.
|
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 834 of file petsc_vector_base.h.
|
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.
|
mutableprivateinherited |
In this map, we count subscriptions for each different identification string supplied to subscribe().
Definition at line 224 of file subscriptor.h.
|
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.
|
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.
|
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.